xref: /petsc/src/snes/impls/fas/fasfunc.c (revision 0dd27c6cfae9dc1af950b720382efbc5314fd943)
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__
365*0dd27c6cSPeter Brune #define __FUNCT__ "SNESFASSetLog"
366*0dd27c6cSPeter Brune /*@
367*0dd27c6cSPeter Brune    SNESFASSetLog - Sets or unsets time logging for various FAS stages on all levels
368*0dd27c6cSPeter Brune 
369*0dd27c6cSPeter Brune    Logically Collective on SNES
370*0dd27c6cSPeter Brune 
371*0dd27c6cSPeter Brune    Input Parameters:
372*0dd27c6cSPeter Brune +  snes   - the FAS context
373*0dd27c6cSPeter Brune -  flg    - monitor or not
374*0dd27c6cSPeter Brune 
375*0dd27c6cSPeter Brune    Level: advanced
376*0dd27c6cSPeter Brune 
377*0dd27c6cSPeter Brune .keywords: FAS, logging
378*0dd27c6cSPeter Brune 
379*0dd27c6cSPeter Brune .seealso: SNESFASSetMonitor()
380*0dd27c6cSPeter Brune @*/
381*0dd27c6cSPeter Brune PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg)
382*0dd27c6cSPeter Brune {
383*0dd27c6cSPeter Brune   SNES_FAS       *fas = (SNES_FAS *)snes->data;
384*0dd27c6cSPeter Brune   PetscErrorCode ierr;
385*0dd27c6cSPeter Brune   PetscBool      isFine;
386*0dd27c6cSPeter Brune   PetscInt       i, levels = fas->levels;
387*0dd27c6cSPeter Brune   SNES           levelsnes;
388*0dd27c6cSPeter Brune   char           eventname[128];
389*0dd27c6cSPeter Brune 
390*0dd27c6cSPeter Brune   PetscFunctionBegin;
391*0dd27c6cSPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
392*0dd27c6cSPeter Brune   if (isFine) {
393*0dd27c6cSPeter Brune     for (i = 0; i < levels; i++) {
394*0dd27c6cSPeter Brune       ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
395*0dd27c6cSPeter Brune       fas = (SNES_FAS *)levelsnes->data;
396*0dd27c6cSPeter Brune       if (flg) {
397*0dd27c6cSPeter Brune         sprintf(eventname,"FASSetup %d",(int)i);
398*0dd27c6cSPeter Brune         ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsetup);CHKERRQ(ierr);
399*0dd27c6cSPeter Brune         sprintf(eventname,"FASSmooth %d",(int)i);
400*0dd27c6cSPeter Brune         ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsolve);CHKERRQ(ierr);
401*0dd27c6cSPeter Brune         sprintf(eventname,"FASResid %d",(int)i);
402*0dd27c6cSPeter Brune         ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventresidual);CHKERRQ(ierr);
403*0dd27c6cSPeter Brune         if (i) {
404*0dd27c6cSPeter Brune           sprintf(eventname,"FASInterp %d",(int)i);
405*0dd27c6cSPeter Brune           ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventinterprestrict);CHKERRQ(ierr);
406*0dd27c6cSPeter Brune         }
407*0dd27c6cSPeter Brune       } else {
408*0dd27c6cSPeter Brune         fas->eventsmoothsetup    = PETSC_FALSE;
409*0dd27c6cSPeter Brune         fas->eventsmoothsolve    = PETSC_FALSE;
410*0dd27c6cSPeter Brune         fas->eventresidual       = PETSC_FALSE;
411*0dd27c6cSPeter Brune         fas->eventinterprestrict = PETSC_FALSE;
412*0dd27c6cSPeter Brune       }
413*0dd27c6cSPeter Brune     }
414*0dd27c6cSPeter Brune   }
415*0dd27c6cSPeter Brune   PetscFunctionReturn(0);
416*0dd27c6cSPeter Brune }
417*0dd27c6cSPeter Brune 
418*0dd27c6cSPeter Brune #undef __FUNCT__
419ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleCreateSmoother_Private"
420ab8d36c9SPeter Brune /*
421ab8d36c9SPeter Brune Creates the default smoother type.
422ab8d36c9SPeter Brune 
42304d7464bSBarry Smith This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level.
424ab8d36c9SPeter Brune 
425ab8d36c9SPeter Brune  */
42622d28d08SBarry Smith PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth)
42722d28d08SBarry Smith {
428ab8d36c9SPeter Brune   SNES_FAS       *fas;
429ab8d36c9SPeter Brune   const char     *optionsprefix;
430ab8d36c9SPeter Brune   char           tprefix[128];
431ab8d36c9SPeter Brune   PetscErrorCode ierr;
432ab8d36c9SPeter Brune   SNES           nsmooth;
43322d28d08SBarry Smith 
434ab8d36c9SPeter Brune   PetscFunctionBegin;
435ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
436ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
437ab8d36c9SPeter Brune   ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr);
438ab8d36c9SPeter Brune   /* create the default smoother */
439ab8d36c9SPeter Brune   ierr = SNESCreate(((PetscObject)snes)->comm, &nsmooth);CHKERRQ(ierr);
440ab8d36c9SPeter Brune   if (fas->level == 0) {
441ab8d36c9SPeter Brune     sprintf(tprefix,"fas_coarse_");
442ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr);
443ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr);
44404d7464bSBarry Smith     ierr = SNESSetType(nsmooth, SNESNEWTONLS);CHKERRQ(ierr);
445e70c42e5SPeter Brune     ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);CHKERRQ(ierr);
446ab8d36c9SPeter Brune   } else {
447ab8d36c9SPeter Brune     sprintf(tprefix,"fas_levels_%d_",(int)fas->level);
448ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr);
449ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr);
450ab8d36c9SPeter Brune     ierr = SNESSetType(nsmooth, SNESNRICHARDSON);CHKERRQ(ierr);
451e70c42e5SPeter Brune     ierr = SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);CHKERRQ(ierr);
452ab8d36c9SPeter Brune   }
453ab8d36c9SPeter Brune   ierr = PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
454f89ba88eSPeter Brune   ierr = PetscLogObjectParent(snes,nsmooth);CHKERRQ(ierr);
455f89ba88eSPeter Brune   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth);CHKERRQ(ierr);
456ab8d36c9SPeter Brune   *smooth = nsmooth;
457ab8d36c9SPeter Brune   PetscFunctionReturn(0);
458ab8d36c9SPeter Brune }
459ab8d36c9SPeter Brune 
460ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */
461ab8d36c9SPeter Brune 
462ab8d36c9SPeter Brune #undef __FUNCT__
463ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleSetCycles"
464ab8d36c9SPeter Brune /*@
465ab8d36c9SPeter Brune    SNESFASCycleSetCycles - Sets the number of cycles on a particular level.
466ab8d36c9SPeter Brune 
467ab8d36c9SPeter Brune    Logically Collective on SNES
468ab8d36c9SPeter Brune 
469ab8d36c9SPeter Brune    Input Parameters:
470ab8d36c9SPeter Brune +  snes   - the multigrid context
471ab8d36c9SPeter Brune .  level  - the level to set the number of cycles on
472ab8d36c9SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
473ab8d36c9SPeter Brune 
474ab8d36c9SPeter Brune    Level: advanced
475ab8d36c9SPeter Brune 
476ab8d36c9SPeter Brune .keywords: SNES, FAS, set, cycles, V-cycle, W-cycle, multigrid
477ab8d36c9SPeter Brune 
478ab8d36c9SPeter Brune .seealso: SNESFASSetCycles()
479ab8d36c9SPeter Brune @*/
48022d28d08SBarry Smith PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles)
48122d28d08SBarry Smith {
482ab8d36c9SPeter Brune   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
483ab8d36c9SPeter Brune   PetscErrorCode ierr;
48422d28d08SBarry Smith 
485ab8d36c9SPeter Brune   PetscFunctionBegin;
486ab8d36c9SPeter Brune   fas->n_cycles = cycles;
487ab8d36c9SPeter Brune   ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr);
488ab8d36c9SPeter Brune   PetscFunctionReturn(0);
489ab8d36c9SPeter Brune }
490ab8d36c9SPeter Brune 
491ab8d36c9SPeter Brune 
492ab8d36c9SPeter Brune #undef __FUNCT__
493ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmoother"
494ab8d36c9SPeter Brune /*@
495ab8d36c9SPeter Brune    SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.
496ab8d36c9SPeter Brune 
497ab8d36c9SPeter Brune    Logically Collective on SNES
498ab8d36c9SPeter Brune 
499ab8d36c9SPeter Brune    Input Parameters:
500ab8d36c9SPeter Brune .  snes   - the multigrid context
501ab8d36c9SPeter Brune 
502ab8d36c9SPeter Brune    Output Parameters:
503ab8d36c9SPeter Brune .  smooth - the smoother
504ab8d36c9SPeter Brune 
505ab8d36c9SPeter Brune    Level: advanced
506ab8d36c9SPeter Brune 
507ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
508ab8d36c9SPeter Brune 
509ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown()
510ab8d36c9SPeter Brune @*/
511ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth)
512ab8d36c9SPeter Brune {
513ab8d36c9SPeter Brune   SNES_FAS       *fas;
514ab8d36c9SPeter Brune   PetscFunctionBegin;
515ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
516ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
517ab8d36c9SPeter Brune   *smooth = fas->smoothd;
518ab8d36c9SPeter Brune   PetscFunctionReturn(0);
519ab8d36c9SPeter Brune }
520ab8d36c9SPeter Brune #undef __FUNCT__
521ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmootherUp"
522ab8d36c9SPeter Brune /*@
523ab8d36c9SPeter Brune    SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level.
524ab8d36c9SPeter Brune 
525ab8d36c9SPeter Brune    Logically Collective on SNES
526ab8d36c9SPeter Brune 
527ab8d36c9SPeter Brune    Input Parameters:
528ab8d36c9SPeter Brune .  snes   - the multigrid context
529ab8d36c9SPeter Brune 
530ab8d36c9SPeter Brune    Output Parameters:
531ab8d36c9SPeter Brune .  smoothu - the smoother
532ab8d36c9SPeter Brune 
533ab8d36c9SPeter Brune    Notes:
534ab8d36c9SPeter Brune    Returns the downsmoother if no up smoother is available.  This enables transparent
535ab8d36c9SPeter Brune    default behavior in the process of the solve.
536ab8d36c9SPeter Brune 
537ab8d36c9SPeter Brune    Level: advanced
538ab8d36c9SPeter Brune 
539ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
540ab8d36c9SPeter Brune 
541ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown()
542ab8d36c9SPeter Brune @*/
543ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu)
544ab8d36c9SPeter Brune {
545ab8d36c9SPeter Brune   SNES_FAS       *fas;
546ab8d36c9SPeter Brune   PetscFunctionBegin;
547ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
548ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
549ab8d36c9SPeter Brune   if (!fas->smoothu) {
550ab8d36c9SPeter Brune     *smoothu = fas->smoothd;
551ab8d36c9SPeter Brune   } else {
552ab8d36c9SPeter Brune     *smoothu = fas->smoothu;
553ab8d36c9SPeter Brune   }
554ab8d36c9SPeter Brune   PetscFunctionReturn(0);
555ab8d36c9SPeter Brune }
556ab8d36c9SPeter Brune 
557ab8d36c9SPeter Brune #undef __FUNCT__
558ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmootherDown"
559ab8d36c9SPeter Brune /*@
560ab8d36c9SPeter Brune    SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level.
561ab8d36c9SPeter Brune 
562ab8d36c9SPeter Brune    Logically Collective on SNES
563ab8d36c9SPeter Brune 
564ab8d36c9SPeter Brune    Input Parameters:
565ab8d36c9SPeter Brune .  snes   - the multigrid context
566ab8d36c9SPeter Brune 
567ab8d36c9SPeter Brune    Output Parameters:
568ab8d36c9SPeter Brune .  smoothd - the smoother
569ab8d36c9SPeter Brune 
570ab8d36c9SPeter Brune    Level: advanced
571ab8d36c9SPeter Brune 
572ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
573ab8d36c9SPeter Brune 
574ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
575ab8d36c9SPeter Brune @*/
576ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd)
577ab8d36c9SPeter Brune {
578ab8d36c9SPeter Brune   SNES_FAS       *fas;
579ab8d36c9SPeter Brune   PetscFunctionBegin;
580ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
581ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
582ab8d36c9SPeter Brune   *smoothd = fas->smoothd;
583ab8d36c9SPeter Brune   PetscFunctionReturn(0);
584ab8d36c9SPeter Brune }
585ab8d36c9SPeter Brune 
586ab8d36c9SPeter Brune 
587ab8d36c9SPeter Brune #undef __FUNCT__
588ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetCorrection"
589ab8d36c9SPeter Brune /*@
590ab8d36c9SPeter Brune    SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level
591ab8d36c9SPeter Brune 
592ab8d36c9SPeter Brune    Logically Collective on SNES
593ab8d36c9SPeter Brune 
594ab8d36c9SPeter Brune    Input Parameters:
595ab8d36c9SPeter Brune .  snes   - the multigrid context
596ab8d36c9SPeter Brune 
597ab8d36c9SPeter Brune    Output Parameters:
598ab8d36c9SPeter Brune .  correction - the coarse correction on this level
599ab8d36c9SPeter Brune 
600ab8d36c9SPeter Brune    Notes:
601ab8d36c9SPeter Brune    Returns PETSC_NULL on the coarsest level.
602ab8d36c9SPeter Brune 
603ab8d36c9SPeter Brune    Level: advanced
604ab8d36c9SPeter Brune 
605ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
606ab8d36c9SPeter Brune 
607ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
608ab8d36c9SPeter Brune @*/
609ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction)
610ab8d36c9SPeter Brune {
611ab8d36c9SPeter Brune   SNES_FAS       *fas;
612ab8d36c9SPeter Brune   PetscFunctionBegin;
613ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
614ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
615ab8d36c9SPeter Brune   *correction = fas->next;
616ab8d36c9SPeter Brune   PetscFunctionReturn(0);
617ab8d36c9SPeter Brune }
618ab8d36c9SPeter Brune 
619ab8d36c9SPeter Brune #undef __FUNCT__
620ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetInterpolation"
621ab8d36c9SPeter Brune /*@
622ab8d36c9SPeter Brune    SNESFASCycleGetInterpolation - Gets the interpolation on this level
623ab8d36c9SPeter Brune 
624ab8d36c9SPeter Brune    Logically Collective on SNES
625ab8d36c9SPeter Brune 
626ab8d36c9SPeter Brune    Input Parameters:
627ab8d36c9SPeter Brune .  snes   - the multigrid context
628ab8d36c9SPeter Brune 
629ab8d36c9SPeter Brune    Output Parameters:
630ab8d36c9SPeter Brune .  mat    - the interpolation operator on this level
631ab8d36c9SPeter Brune 
632ab8d36c9SPeter Brune    Level: developer
633ab8d36c9SPeter Brune 
634ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
635ab8d36c9SPeter Brune 
636ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
637ab8d36c9SPeter Brune @*/
638ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat)
639ab8d36c9SPeter Brune {
640ab8d36c9SPeter Brune   SNES_FAS       *fas;
641ab8d36c9SPeter Brune   PetscFunctionBegin;
642ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
643ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
644ab8d36c9SPeter Brune   *mat = fas->interpolate;
645ab8d36c9SPeter Brune   PetscFunctionReturn(0);
646ab8d36c9SPeter Brune }
647ab8d36c9SPeter Brune 
648ab8d36c9SPeter Brune 
649ab8d36c9SPeter Brune #undef __FUNCT__
650ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetRestriction"
651ab8d36c9SPeter Brune /*@
652ab8d36c9SPeter Brune    SNESFASCycleGetRestriction - Gets the restriction on this level
653ab8d36c9SPeter Brune 
654ab8d36c9SPeter Brune    Logically Collective on SNES
655ab8d36c9SPeter Brune 
656ab8d36c9SPeter Brune    Input Parameters:
657ab8d36c9SPeter Brune .  snes   - the multigrid context
658ab8d36c9SPeter Brune 
659ab8d36c9SPeter Brune    Output Parameters:
660ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
661ab8d36c9SPeter Brune 
662ab8d36c9SPeter Brune    Level: developer
663ab8d36c9SPeter Brune 
664ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
665ab8d36c9SPeter Brune 
666ab8d36c9SPeter Brune .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation()
667ab8d36c9SPeter Brune @*/
668ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat)
669ab8d36c9SPeter Brune {
670ab8d36c9SPeter Brune   SNES_FAS       *fas;
671ab8d36c9SPeter Brune   PetscFunctionBegin;
672ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
673ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
674ab8d36c9SPeter Brune   *mat = fas->restrct;
675ab8d36c9SPeter Brune   PetscFunctionReturn(0);
676ab8d36c9SPeter Brune }
677ab8d36c9SPeter Brune 
678ab8d36c9SPeter Brune 
679ab8d36c9SPeter Brune #undef __FUNCT__
680ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetInjection"
681ab8d36c9SPeter Brune /*@
682ab8d36c9SPeter Brune    SNESFASCycleGetInjection - Gets the injection on this level
683ab8d36c9SPeter Brune 
684ab8d36c9SPeter Brune    Logically Collective on SNES
685ab8d36c9SPeter Brune 
686ab8d36c9SPeter Brune    Input Parameters:
687ab8d36c9SPeter Brune .  snes   - the multigrid context
688ab8d36c9SPeter Brune 
689ab8d36c9SPeter Brune    Output Parameters:
690ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
691ab8d36c9SPeter Brune 
692ab8d36c9SPeter Brune    Level: developer
693ab8d36c9SPeter Brune 
694ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
695ab8d36c9SPeter Brune 
696ab8d36c9SPeter Brune .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction()
697ab8d36c9SPeter Brune @*/
698ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat)
699ab8d36c9SPeter Brune {
700ab8d36c9SPeter Brune   SNES_FAS       *fas;
701ab8d36c9SPeter Brune   PetscFunctionBegin;
702ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
703ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
704ab8d36c9SPeter Brune   *mat = fas->inject;
705ab8d36c9SPeter Brune   PetscFunctionReturn(0);
706ab8d36c9SPeter Brune }
707ab8d36c9SPeter Brune 
708ab8d36c9SPeter Brune #undef __FUNCT__
709ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetRScale"
710ab8d36c9SPeter Brune /*@
711ab8d36c9SPeter Brune    SNESFASCycleGetRScale - Gets the injection on this level
712ab8d36c9SPeter Brune 
713ab8d36c9SPeter Brune    Logically Collective on SNES
714ab8d36c9SPeter Brune 
715ab8d36c9SPeter Brune    Input Parameters:
716ab8d36c9SPeter Brune .  snes   - the multigrid context
717ab8d36c9SPeter Brune 
718ab8d36c9SPeter Brune    Output Parameters:
719ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
720ab8d36c9SPeter Brune 
721ab8d36c9SPeter Brune    Level: developer
722ab8d36c9SPeter Brune 
723ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
724ab8d36c9SPeter Brune 
725ab8d36c9SPeter Brune .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale()
726ab8d36c9SPeter Brune @*/
727ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec)
728ab8d36c9SPeter Brune {
729ab8d36c9SPeter Brune   SNES_FAS       *fas;
730ab8d36c9SPeter Brune   PetscFunctionBegin;
731ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
732ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
733ab8d36c9SPeter Brune   *vec = fas->rscale;
734ab8d36c9SPeter Brune   PetscFunctionReturn(0);
735ab8d36c9SPeter Brune }
736ab8d36c9SPeter Brune 
737ab8d36c9SPeter Brune #undef __FUNCT__
738ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleIsFine"
739ab8d36c9SPeter Brune /*@
740ab8d36c9SPeter Brune    SNESFASCycleIsFine - Determines if a given cycle is the fine level.
741ab8d36c9SPeter Brune 
742ab8d36c9SPeter Brune    Logically Collective on SNES
743ab8d36c9SPeter Brune 
744ab8d36c9SPeter Brune    Input Parameters:
745ab8d36c9SPeter Brune .  snes   - the FAS context
746ab8d36c9SPeter Brune 
747ab8d36c9SPeter Brune    Output Parameters:
748ab8d36c9SPeter Brune .  flg - indicates if this is the fine level or not
749ab8d36c9SPeter Brune 
750ab8d36c9SPeter Brune    Level: advanced
751ab8d36c9SPeter Brune 
752ab8d36c9SPeter Brune .keywords: SNES, FAS
753ab8d36c9SPeter Brune 
754ab8d36c9SPeter Brune .seealso: SNESFASSetLevels()
755ab8d36c9SPeter Brune @*/
756ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg)
757ab8d36c9SPeter Brune {
758ab8d36c9SPeter Brune   SNES_FAS       *fas;
759ab8d36c9SPeter Brune   PetscFunctionBegin;
760ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
761ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
762ab8d36c9SPeter Brune   if (fas->level == fas->levels - 1) {
763ab8d36c9SPeter Brune     *flg = PETSC_TRUE;
764ab8d36c9SPeter Brune   } else {
765ab8d36c9SPeter Brune     *flg = PETSC_FALSE;
766ab8d36c9SPeter Brune   }
767ab8d36c9SPeter Brune   PetscFunctionReturn(0);
768ab8d36c9SPeter Brune }
769ab8d36c9SPeter Brune 
770ab8d36c9SPeter Brune /* ---------- functions called on the finest level that return level-specific information ---------- */
771ab8d36c9SPeter Brune 
772ab8d36c9SPeter Brune #undef __FUNCT__
773ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
774ab8d36c9SPeter Brune /*@
775ab8d36c9SPeter Brune    SNESFASSetInterpolation - Sets the function to be used to calculate the
776ab8d36c9SPeter Brune    interpolation from l-1 to the lth level
777ab8d36c9SPeter Brune 
778ab8d36c9SPeter Brune    Input Parameters:
779ab8d36c9SPeter Brune +  snes      - the multigrid context
780ab8d36c9SPeter Brune .  mat       - the interpolation operator
781ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
782ab8d36c9SPeter Brune 
783ab8d36c9SPeter Brune    Level: advanced
784ab8d36c9SPeter Brune 
785ab8d36c9SPeter Brune    Notes:
786ab8d36c9SPeter Brune           Usually this is the same matrix used also to set the restriction
787ab8d36c9SPeter Brune     for the same level.
788ab8d36c9SPeter Brune 
789ab8d36c9SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
790ab8d36c9SPeter Brune     out from the matrix size which one it is.
791ab8d36c9SPeter Brune 
792ab8d36c9SPeter Brune .keywords:  FAS, multigrid, set, interpolate, level
793ab8d36c9SPeter Brune 
794ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
795ab8d36c9SPeter Brune @*/
7960adebc6cSBarry Smith PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat)
7970adebc6cSBarry Smith {
79822d28d08SBarry Smith   SNES_FAS       *fas;
799ab8d36c9SPeter Brune   PetscErrorCode ierr;
800ab8d36c9SPeter Brune   SNES           levelsnes;
80122d28d08SBarry Smith 
802ab8d36c9SPeter Brune   PetscFunctionBegin;
803ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
804ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
805ab8d36c9SPeter Brune   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
806ab8d36c9SPeter Brune   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
807ab8d36c9SPeter Brune   fas->interpolate = mat;
808ab8d36c9SPeter Brune   PetscFunctionReturn(0);
809ab8d36c9SPeter Brune }
810ab8d36c9SPeter Brune 
811ab8d36c9SPeter Brune #undef __FUNCT__
812ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetInterpolation"
813ab8d36c9SPeter Brune /*@
814ab8d36c9SPeter Brune    SNESFASGetInterpolation - Gets the matrix used to calculate the
815ab8d36c9SPeter Brune    interpolation from l-1 to the lth level
816ab8d36c9SPeter Brune 
817ab8d36c9SPeter Brune    Input Parameters:
818ab8d36c9SPeter Brune +  snes      - the multigrid context
819ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
820ab8d36c9SPeter Brune 
821ab8d36c9SPeter Brune    Output Parameters:
822ab8d36c9SPeter Brune .  mat       - the interpolation operator
823ab8d36c9SPeter Brune 
824ab8d36c9SPeter Brune    Level: advanced
825ab8d36c9SPeter Brune 
826ab8d36c9SPeter Brune .keywords:  FAS, multigrid, get, interpolate, level
827ab8d36c9SPeter Brune 
828ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale()
829ab8d36c9SPeter Brune @*/
83022d28d08SBarry Smith PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat)
83122d28d08SBarry Smith {
83222d28d08SBarry Smith   SNES_FAS       *fas;
833ab8d36c9SPeter Brune   PetscErrorCode ierr;
834ab8d36c9SPeter Brune   SNES           levelsnes;
83522d28d08SBarry Smith 
836ab8d36c9SPeter Brune   PetscFunctionBegin;
837ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
838ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
839ab8d36c9SPeter Brune   *mat = fas->interpolate;
840ab8d36c9SPeter Brune   PetscFunctionReturn(0);
841ab8d36c9SPeter Brune }
842ab8d36c9SPeter Brune 
843ab8d36c9SPeter Brune #undef __FUNCT__
844ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
845ab8d36c9SPeter Brune /*@
846ab8d36c9SPeter Brune    SNESFASSetRestriction - Sets the function to be used to restrict the defect
847ab8d36c9SPeter Brune    from level l to l-1.
848ab8d36c9SPeter Brune 
849ab8d36c9SPeter Brune    Input Parameters:
850ab8d36c9SPeter Brune +  snes  - the multigrid context
851ab8d36c9SPeter Brune .  mat   - the restriction matrix
852ab8d36c9SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
853ab8d36c9SPeter Brune 
854ab8d36c9SPeter Brune    Level: advanced
855ab8d36c9SPeter Brune 
856ab8d36c9SPeter Brune    Notes:
857ab8d36c9SPeter Brune           Usually this is the same matrix used also to set the interpolation
858ab8d36c9SPeter Brune     for the same level.
859ab8d36c9SPeter Brune 
860ab8d36c9SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
861ab8d36c9SPeter Brune     out from the matrix size which one it is.
862ab8d36c9SPeter Brune 
863ab8d36c9SPeter Brune          If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
864ab8d36c9SPeter Brune     is used.
865ab8d36c9SPeter Brune 
866ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
867ab8d36c9SPeter Brune 
868ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
869ab8d36c9SPeter Brune @*/
87022d28d08SBarry Smith PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat)
87122d28d08SBarry Smith {
87222d28d08SBarry Smith   SNES_FAS       *fas;
873ab8d36c9SPeter Brune   PetscErrorCode ierr;
874ab8d36c9SPeter Brune   SNES           levelsnes;
87522d28d08SBarry Smith 
876ab8d36c9SPeter Brune   PetscFunctionBegin;
877ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
878ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
879ab8d36c9SPeter Brune   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
880ab8d36c9SPeter Brune   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
881ab8d36c9SPeter Brune   fas->restrct = mat;
882ab8d36c9SPeter Brune   PetscFunctionReturn(0);
883ab8d36c9SPeter Brune }
884ab8d36c9SPeter Brune 
885ab8d36c9SPeter Brune #undef __FUNCT__
886ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetRestriction"
887ab8d36c9SPeter Brune /*@
888ab8d36c9SPeter Brune    SNESFASGetRestriction - Gets the matrix used to calculate the
889ab8d36c9SPeter Brune    restriction from l to the l-1th level
890ab8d36c9SPeter Brune 
891ab8d36c9SPeter Brune    Input Parameters:
892ab8d36c9SPeter Brune +  snes      - the multigrid context
893ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
894ab8d36c9SPeter Brune 
895ab8d36c9SPeter Brune    Output Parameters:
896ab8d36c9SPeter Brune .  mat       - the interpolation operator
897ab8d36c9SPeter Brune 
898ab8d36c9SPeter Brune    Level: advanced
899ab8d36c9SPeter Brune 
900ab8d36c9SPeter Brune .keywords:  FAS, multigrid, get, restrict, level
901ab8d36c9SPeter Brune 
902ab8d36c9SPeter Brune .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale()
903ab8d36c9SPeter Brune @*/
90422d28d08SBarry Smith PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat)
90522d28d08SBarry Smith {
90622d28d08SBarry Smith   SNES_FAS       *fas;
907ab8d36c9SPeter Brune   PetscErrorCode ierr;
908ab8d36c9SPeter Brune   SNES           levelsnes;
90922d28d08SBarry Smith 
910ab8d36c9SPeter Brune   PetscFunctionBegin;
911ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
912ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
913ab8d36c9SPeter Brune   *mat = fas->restrct;
914ab8d36c9SPeter Brune   PetscFunctionReturn(0);
915ab8d36c9SPeter Brune }
916ab8d36c9SPeter Brune 
917ab8d36c9SPeter Brune 
918ab8d36c9SPeter Brune #undef __FUNCT__
919ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInjection"
920ab8d36c9SPeter Brune /*@
921ab8d36c9SPeter Brune    SNESFASSetInjection - Sets the function to be used to inject the solution
922ab8d36c9SPeter Brune    from level l to l-1.
923ab8d36c9SPeter Brune 
924ab8d36c9SPeter Brune    Input Parameters:
925ab8d36c9SPeter Brune  +  snes  - the multigrid context
926ab8d36c9SPeter Brune .  mat   - the restriction matrix
927ab8d36c9SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
928ab8d36c9SPeter Brune 
929ab8d36c9SPeter Brune    Level: advanced
930ab8d36c9SPeter Brune 
931ab8d36c9SPeter Brune    Notes:
932ab8d36c9SPeter Brune          If you do not set this, the restriction and rscale is used to
933ab8d36c9SPeter Brune    project the solution instead.
934ab8d36c9SPeter Brune 
935ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
936ab8d36c9SPeter Brune 
937ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
938ab8d36c9SPeter Brune @*/
93922d28d08SBarry Smith PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat)
94022d28d08SBarry Smith {
94122d28d08SBarry Smith   SNES_FAS       *fas;
942ab8d36c9SPeter Brune   PetscErrorCode ierr;
943ab8d36c9SPeter Brune   SNES           levelsnes;
94422d28d08SBarry Smith 
945ab8d36c9SPeter Brune   PetscFunctionBegin;
946ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
947ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
948ab8d36c9SPeter Brune   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
949ab8d36c9SPeter Brune   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
950ab8d36c9SPeter Brune   fas->inject = mat;
951ab8d36c9SPeter Brune   PetscFunctionReturn(0);
952ab8d36c9SPeter Brune }
953ab8d36c9SPeter Brune 
954ab8d36c9SPeter Brune 
955ab8d36c9SPeter Brune #undef __FUNCT__
956ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetInjection"
957ab8d36c9SPeter Brune /*@
958ab8d36c9SPeter Brune    SNESFASGetInjection - Gets the matrix used to calculate the
959ab8d36c9SPeter Brune    injection from l-1 to the lth level
960ab8d36c9SPeter Brune 
961ab8d36c9SPeter Brune    Input Parameters:
962ab8d36c9SPeter Brune +  snes      - the multigrid context
963ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
964ab8d36c9SPeter Brune 
965ab8d36c9SPeter Brune    Output Parameters:
966ab8d36c9SPeter Brune .  mat       - the injection operator
967ab8d36c9SPeter Brune 
968ab8d36c9SPeter Brune    Level: advanced
969ab8d36c9SPeter Brune 
970ab8d36c9SPeter Brune .keywords:  FAS, multigrid, get, restrict, level
971ab8d36c9SPeter Brune 
972ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale()
973ab8d36c9SPeter Brune @*/
97422d28d08SBarry Smith PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat)
97522d28d08SBarry Smith {
97622d28d08SBarry Smith   SNES_FAS       *fas;
977ab8d36c9SPeter Brune   PetscErrorCode ierr;
978ab8d36c9SPeter Brune   SNES           levelsnes;
97922d28d08SBarry Smith 
980ab8d36c9SPeter Brune   PetscFunctionBegin;
981ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
982ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
983ab8d36c9SPeter Brune   *mat = fas->inject;
984ab8d36c9SPeter Brune   PetscFunctionReturn(0);
985ab8d36c9SPeter Brune }
986ab8d36c9SPeter Brune 
987ab8d36c9SPeter Brune #undef __FUNCT__
988ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
989ab8d36c9SPeter Brune /*@
990ab8d36c9SPeter Brune    SNESFASSetRScale - Sets the scaling factor of the restriction
991ab8d36c9SPeter Brune    operator from level l to l-1.
992ab8d36c9SPeter Brune 
993ab8d36c9SPeter Brune    Input Parameters:
994ab8d36c9SPeter Brune +  snes   - the multigrid context
995ab8d36c9SPeter Brune .  rscale - the restriction scaling
996ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply [Do not supply 0]
997ab8d36c9SPeter Brune 
998ab8d36c9SPeter Brune    Level: advanced
999ab8d36c9SPeter Brune 
1000ab8d36c9SPeter Brune    Notes:
1001ab8d36c9SPeter Brune          This is only used in the case that the injection is not set.
1002ab8d36c9SPeter Brune 
1003ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
1004ab8d36c9SPeter Brune 
1005ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1006ab8d36c9SPeter Brune @*/
100722d28d08SBarry Smith PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale)
100822d28d08SBarry Smith {
1009ab8d36c9SPeter Brune   SNES_FAS       *fas;
1010ab8d36c9SPeter Brune   PetscErrorCode ierr;
1011ab8d36c9SPeter Brune   SNES           levelsnes;
101222d28d08SBarry Smith 
1013ab8d36c9SPeter Brune   PetscFunctionBegin;
1014ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1015ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
1016ab8d36c9SPeter Brune   ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr);
1017ab8d36c9SPeter Brune   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
1018ab8d36c9SPeter Brune   fas->rscale = rscale;
1019ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1020ab8d36c9SPeter Brune }
1021ab8d36c9SPeter Brune 
1022ab8d36c9SPeter Brune #undef __FUNCT__
1023ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmoother"
1024ab8d36c9SPeter Brune /*@
1025ab8d36c9SPeter Brune    SNESFASGetSmoother - Gets the default smoother on a level.
1026ab8d36c9SPeter Brune 
1027ab8d36c9SPeter Brune    Input Parameters:
1028ab8d36c9SPeter Brune +  snes   - the multigrid context
1029ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply
1030ab8d36c9SPeter Brune 
1031ab8d36c9SPeter Brune    Output Parameters:
1032ab8d36c9SPeter Brune    smooth  - the smoother
1033ab8d36c9SPeter Brune 
1034ab8d36c9SPeter Brune    Level: advanced
1035ab8d36c9SPeter Brune 
1036ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level
1037ab8d36c9SPeter Brune 
1038ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1039ab8d36c9SPeter Brune @*/
104022d28d08SBarry Smith PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth)
104122d28d08SBarry Smith {
1042ab8d36c9SPeter Brune   SNES_FAS       *fas;
1043ab8d36c9SPeter Brune   PetscErrorCode ierr;
1044ab8d36c9SPeter Brune   SNES           levelsnes;
104522d28d08SBarry Smith 
1046ab8d36c9SPeter Brune   PetscFunctionBegin;
1047ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1048ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
1049ab8d36c9SPeter Brune   if (!fas->smoothd) {
1050ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
1051ab8d36c9SPeter Brune   }
1052ab8d36c9SPeter Brune   *smooth = fas->smoothd;
1053ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1054ab8d36c9SPeter Brune }
1055ab8d36c9SPeter Brune 
1056ab8d36c9SPeter Brune #undef __FUNCT__
1057ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherDown"
1058ab8d36c9SPeter Brune /*@
1059ab8d36c9SPeter Brune    SNESFASGetSmootherDown - Gets the downsmoother on a level.
1060ab8d36c9SPeter Brune 
1061ab8d36c9SPeter Brune    Input Parameters:
1062ab8d36c9SPeter Brune +  snes   - the multigrid context
1063ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply
1064ab8d36c9SPeter Brune 
1065ab8d36c9SPeter Brune    Output Parameters:
1066ab8d36c9SPeter Brune    smooth  - the smoother
1067ab8d36c9SPeter Brune 
1068ab8d36c9SPeter Brune    Level: advanced
1069ab8d36c9SPeter Brune 
1070ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level
1071ab8d36c9SPeter Brune 
1072ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1073ab8d36c9SPeter Brune @*/
107422d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth)
107522d28d08SBarry Smith {
1076ab8d36c9SPeter Brune   SNES_FAS       *fas;
1077ab8d36c9SPeter Brune   PetscErrorCode ierr;
1078ab8d36c9SPeter Brune   SNES           levelsnes;
107922d28d08SBarry Smith 
1080ab8d36c9SPeter Brune   PetscFunctionBegin;
1081ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1082ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
1083ab8d36c9SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1084ab8d36c9SPeter Brune   if (!fas->smoothd) {
1085ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
1086ab8d36c9SPeter Brune   }
1087ab8d36c9SPeter Brune   if (!fas->smoothu) {
1088ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr);
1089ab8d36c9SPeter Brune   }
1090ab8d36c9SPeter Brune   *smooth = fas->smoothd;
1091ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1092ab8d36c9SPeter Brune }
1093ab8d36c9SPeter Brune 
1094ab8d36c9SPeter Brune #undef __FUNCT__
1095ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherUp"
1096ab8d36c9SPeter Brune /*@
1097ab8d36c9SPeter Brune    SNESFASGetSmootherUp - Gets the upsmoother on a level.
1098ab8d36c9SPeter Brune 
1099ab8d36c9SPeter Brune    Input Parameters:
1100ab8d36c9SPeter Brune +  snes   - the multigrid context
1101ab8d36c9SPeter Brune -  level  - the level (0 is coarsest)
1102ab8d36c9SPeter Brune 
1103ab8d36c9SPeter Brune    Output Parameters:
1104ab8d36c9SPeter Brune    smooth  - the smoother
1105ab8d36c9SPeter Brune 
1106ab8d36c9SPeter Brune    Level: advanced
1107ab8d36c9SPeter Brune 
1108ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level
1109ab8d36c9SPeter Brune 
1110ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1111ab8d36c9SPeter Brune @*/
111222d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth)
111322d28d08SBarry Smith {
1114ab8d36c9SPeter Brune   SNES_FAS       *fas;
1115ab8d36c9SPeter Brune   PetscErrorCode ierr;
1116ab8d36c9SPeter Brune   SNES           levelsnes;
111722d28d08SBarry Smith 
1118ab8d36c9SPeter Brune   PetscFunctionBegin;
1119ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1120ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
1121ab8d36c9SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1122ab8d36c9SPeter Brune   if (!fas->smoothd) {
1123ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
1124ab8d36c9SPeter Brune   }
1125ab8d36c9SPeter Brune   if (!fas->smoothu) {
1126ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr);
1127ab8d36c9SPeter Brune   }
1128ab8d36c9SPeter Brune   *smooth = fas->smoothu;
1129ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1130ab8d36c9SPeter Brune }
1131d6ad1212SPeter Brune 
1132d6ad1212SPeter Brune #undef __FUNCT__
1133d6ad1212SPeter Brune #define __FUNCT__ "SNESFASGetCoarseSolve"
1134d6ad1212SPeter Brune /*@
1135d6ad1212SPeter Brune    SNESFASGetCoarseSolve - Gets the coarsest solver.
1136d6ad1212SPeter Brune 
1137d6ad1212SPeter Brune    Input Parameters:
1138d6ad1212SPeter Brune +  snes   - the multigrid context
1139d6ad1212SPeter Brune 
1140d6ad1212SPeter Brune    Output Parameters:
1141d6ad1212SPeter Brune    solve  - the coarse-level solver
1142d6ad1212SPeter Brune 
1143d6ad1212SPeter Brune    Level: advanced
1144d6ad1212SPeter Brune 
1145d6ad1212SPeter Brune .keywords: FAS, MG, get, multigrid, solver, coarse
1146d6ad1212SPeter Brune 
1147d6ad1212SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1148d6ad1212SPeter Brune @*/
114922d28d08SBarry Smith PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *smooth)
115022d28d08SBarry Smith {
1151d6ad1212SPeter Brune   SNES_FAS       *fas;
1152d6ad1212SPeter Brune   PetscErrorCode ierr;
1153d6ad1212SPeter Brune   SNES           levelsnes;
115422d28d08SBarry Smith 
1155d6ad1212SPeter Brune   PetscFunctionBegin;
1156d6ad1212SPeter Brune   ierr = SNESFASGetCycleSNES(snes, 0, &levelsnes);CHKERRQ(ierr);
1157d6ad1212SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
1158d6ad1212SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1159d6ad1212SPeter Brune   if (!fas->smoothd) {
1160d6ad1212SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
1161d6ad1212SPeter Brune   }
1162d6ad1212SPeter Brune   *smooth = fas->smoothd;
1163d6ad1212SPeter Brune   PetscFunctionReturn(0);
1164d6ad1212SPeter Brune }
1165