xref: /petsc/src/snes/impls/fas/fasfunc.c (revision 3bb1ff401821b9e2ae019d3e61bc8ab4bd4e59d5)
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
770298fd71SBarry Smith             problems on smaller sets of processors. Use 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;
101ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
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);
1090298fd71SBarry Smith   fas->next     = NULL;
1100298fd71SBarry Smith   fas->previous = 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;
1190298fd71SBarry Smith     fas->next   = 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);
1291aa26658SKarl Rupp 
130ab8d36c9SPeter Brune       ((SNES_FAS*)fas->next->data)->previous = prevsnes;
1311aa26658SKarl Rupp 
132ab8d36c9SPeter Brune       prevsnes = fas->next;
133ab8d36c9SPeter Brune       fas      = (SNES_FAS*)prevsnes->data;
134ab8d36c9SPeter Brune     }
135ab8d36c9SPeter Brune   }
136ab8d36c9SPeter Brune   PetscFunctionReturn(0);
137ab8d36c9SPeter Brune }
138ab8d36c9SPeter Brune 
139ab8d36c9SPeter Brune 
140ab8d36c9SPeter Brune #undef __FUNCT__
141ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetLevels"
142ab8d36c9SPeter Brune /*@
143ab8d36c9SPeter Brune    SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids
144ab8d36c9SPeter Brune 
145ab8d36c9SPeter Brune    Input Parameter:
146ab8d36c9SPeter Brune .  snes - the nonlinear solver context
147ab8d36c9SPeter Brune 
148ab8d36c9SPeter Brune    Output parameter:
149ab8d36c9SPeter Brune .  levels - the number of levels
150ab8d36c9SPeter Brune 
151ab8d36c9SPeter Brune    Level: advanced
152ab8d36c9SPeter Brune 
153ab8d36c9SPeter Brune .keywords: MG, get, levels, multigrid
154ab8d36c9SPeter Brune 
155ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels()
156ab8d36c9SPeter Brune @*/
15722d28d08SBarry Smith PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels)
15822d28d08SBarry Smith {
159ab8d36c9SPeter Brune   SNES_FAS * fas = (SNES_FAS*)snes->data;
1605fd66863SKarl Rupp 
161ab8d36c9SPeter Brune   PetscFunctionBegin;
162ab8d36c9SPeter Brune   *levels = fas->levels;
163ab8d36c9SPeter Brune   PetscFunctionReturn(0);
164ab8d36c9SPeter Brune }
165ab8d36c9SPeter Brune 
166ab8d36c9SPeter Brune 
167ab8d36c9SPeter Brune #undef __FUNCT__
168ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetCycleSNES"
169ab8d36c9SPeter Brune /*@
170ab8d36c9SPeter Brune    SNESFASGetCycleSNES - Gets the SNES corresponding to a particular
171ab8d36c9SPeter Brune    level of the FAS hierarchy.
172ab8d36c9SPeter Brune 
173ab8d36c9SPeter Brune    Input Parameters:
174ab8d36c9SPeter Brune +  snes    - the multigrid context
175ab8d36c9SPeter Brune    level   - the level to get
176ab8d36c9SPeter Brune -  lsnes   - whether to use the nonlinear smoother or not
177ab8d36c9SPeter Brune 
178ab8d36c9SPeter Brune    Level: advanced
179ab8d36c9SPeter Brune 
180ab8d36c9SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
181ab8d36c9SPeter Brune 
182ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels()
183ab8d36c9SPeter Brune @*/
18422d28d08SBarry Smith PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes)
18522d28d08SBarry Smith {
186ab8d36c9SPeter Brune   SNES_FAS *fas = (SNES_FAS*)snes->data;
187ab8d36c9SPeter Brune   PetscInt i;
188ab8d36c9SPeter Brune 
189ab8d36c9SPeter Brune   PetscFunctionBegin;
190ce94432eSBarry Smith   if (level > fas->levels-1) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels);
191ce94432eSBarry Smith   if (fas->level !=  fas->levels - 1) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"SNESFASGetCycleSNES may only be called on the finest-level SNES.",level,fas->level);
192ab8d36c9SPeter Brune 
193ab8d36c9SPeter Brune   *lsnes = snes;
194ab8d36c9SPeter Brune   for (i = fas->level; i > level; i--) {
195ab8d36c9SPeter Brune     *lsnes = fas->next;
196ab8d36c9SPeter Brune     fas    = (SNES_FAS*)(*lsnes)->data;
197ab8d36c9SPeter Brune   }
198ce94432eSBarry Smith   if (fas->level != level) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
199ab8d36c9SPeter Brune   PetscFunctionReturn(0);
200ab8d36c9SPeter Brune }
201ab8d36c9SPeter Brune 
202ab8d36c9SPeter Brune #undef __FUNCT__
203ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp"
204ab8d36c9SPeter Brune /*@
205ab8d36c9SPeter Brune    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
206ab8d36c9SPeter Brune    use on all levels.
207ab8d36c9SPeter Brune 
208ab8d36c9SPeter Brune    Logically Collective on SNES
209ab8d36c9SPeter Brune 
210ab8d36c9SPeter Brune    Input Parameters:
211ab8d36c9SPeter Brune +  snes - the multigrid context
212ab8d36c9SPeter Brune -  n    - the number of smoothing steps
213ab8d36c9SPeter Brune 
214ab8d36c9SPeter Brune    Options Database Key:
215ab8d36c9SPeter Brune .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
216ab8d36c9SPeter Brune 
217ab8d36c9SPeter Brune    Level: advanced
218ab8d36c9SPeter Brune 
219ab8d36c9SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
220ab8d36c9SPeter Brune 
221ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothDown()
222ab8d36c9SPeter Brune @*/
22322d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n)
22422d28d08SBarry Smith {
225ab8d36c9SPeter Brune   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
22622d28d08SBarry Smith   PetscErrorCode ierr;
22722d28d08SBarry Smith 
228ab8d36c9SPeter Brune   PetscFunctionBegin;
229ab8d36c9SPeter Brune   fas->max_up_it = n;
230656ede7eSPeter Brune   if (!fas->smoothu && fas->level != 0) {
23122d28d08SBarry Smith     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr);
232ab8d36c9SPeter Brune   }
23322d28d08SBarry Smith   if (fas->smoothu) {
23422d28d08SBarry Smith     ierr = SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs);CHKERRQ(ierr);
23522d28d08SBarry Smith   }
236ab8d36c9SPeter Brune   if (fas->next) {
237ab8d36c9SPeter Brune     ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr);
238ab8d36c9SPeter Brune   }
239ab8d36c9SPeter Brune   PetscFunctionReturn(0);
240ab8d36c9SPeter Brune }
241ab8d36c9SPeter Brune 
242ab8d36c9SPeter Brune #undef __FUNCT__
243ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown"
244ab8d36c9SPeter Brune /*@
245ab8d36c9SPeter Brune    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
246ab8d36c9SPeter Brune    use on all levels.
247ab8d36c9SPeter Brune 
248ab8d36c9SPeter Brune    Logically Collective on SNES
249ab8d36c9SPeter Brune 
250ab8d36c9SPeter Brune    Input Parameters:
251ab8d36c9SPeter Brune +  snes - the multigrid context
252ab8d36c9SPeter Brune -  n    - the number of smoothing steps
253ab8d36c9SPeter Brune 
254ab8d36c9SPeter Brune    Options Database Key:
255ab8d36c9SPeter Brune .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
256ab8d36c9SPeter Brune 
257ab8d36c9SPeter Brune    Level: advanced
258ab8d36c9SPeter Brune 
259ab8d36c9SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
260ab8d36c9SPeter Brune 
261ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothUp()
262ab8d36c9SPeter Brune @*/
26322d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n)
26422d28d08SBarry Smith {
265ab8d36c9SPeter Brune   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
266ab8d36c9SPeter Brune   PetscErrorCode ierr = 0;
26722d28d08SBarry Smith 
268ab8d36c9SPeter Brune   PetscFunctionBegin;
269ab8d36c9SPeter Brune   if (!fas->smoothd) {
27022d28d08SBarry Smith     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
271ab8d36c9SPeter Brune   }
272ab8d36c9SPeter Brune   ierr = SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs);CHKERRQ(ierr);
2731aa26658SKarl Rupp 
274ab8d36c9SPeter Brune   fas->max_down_it = n;
275ab8d36c9SPeter Brune   if (fas->next) {
276ab8d36c9SPeter Brune     ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr);
277ab8d36c9SPeter Brune   }
278ab8d36c9SPeter Brune   PetscFunctionReturn(0);
279ab8d36c9SPeter Brune }
280ab8d36c9SPeter Brune 
281ab8d36c9SPeter Brune 
282ab8d36c9SPeter Brune #undef __FUNCT__
283ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetCycles"
284ab8d36c9SPeter Brune /*@
285ab8d36c9SPeter Brune    SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited.  Use SNESFASSetCyclesOnLevel() for more
286ab8d36c9SPeter Brune    complicated cycling.
287ab8d36c9SPeter Brune 
288ab8d36c9SPeter Brune    Logically Collective on SNES
289ab8d36c9SPeter Brune 
290ab8d36c9SPeter Brune    Input Parameters:
291ab8d36c9SPeter Brune +  snes   - the multigrid context
292ab8d36c9SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
293ab8d36c9SPeter Brune 
294ab8d36c9SPeter Brune    Options Database Key:
295ab8d36c9SPeter Brune $  -snes_fas_cycles 1 or 2
296ab8d36c9SPeter Brune 
297ab8d36c9SPeter Brune    Level: advanced
298ab8d36c9SPeter Brune 
299ab8d36c9SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
300ab8d36c9SPeter Brune 
301ab8d36c9SPeter Brune .seealso: SNESFASSetCyclesOnLevel()
302ab8d36c9SPeter Brune @*/
30322d28d08SBarry Smith PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles)
30422d28d08SBarry Smith {
305ab8d36c9SPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
306ab8d36c9SPeter Brune   PetscErrorCode ierr;
307ab8d36c9SPeter Brune   PetscBool      isFine;
30822d28d08SBarry Smith 
309ab8d36c9SPeter Brune   PetscFunctionBegin;
31022d28d08SBarry Smith   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
3111aa26658SKarl Rupp 
312ab8d36c9SPeter Brune   fas->n_cycles = cycles;
3131aa26658SKarl Rupp   if (!isFine) {
314ab8d36c9SPeter Brune     ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr);
3151aa26658SKarl Rupp   }
316ab8d36c9SPeter Brune   if (fas->next) {
317ab8d36c9SPeter Brune     ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr);
318ab8d36c9SPeter Brune   }
319ab8d36c9SPeter Brune   PetscFunctionReturn(0);
320ab8d36c9SPeter Brune }
321ab8d36c9SPeter Brune 
322c8c899caSPeter Brune 
323c8c899caSPeter Brune #undef __FUNCT__
324c8c899caSPeter Brune #define __FUNCT__ "SNESFASSetMonitor"
325c8c899caSPeter Brune /*@
326c8c899caSPeter Brune    SNESFASSetMonitor - Sets the method-specific cycle monitoring
327c8c899caSPeter Brune 
328c8c899caSPeter Brune    Logically Collective on SNES
329c8c899caSPeter Brune 
330c8c899caSPeter Brune    Input Parameters:
331c8c899caSPeter Brune +  snes   - the FAS context
332c8c899caSPeter Brune -  flg    - monitor or not
333c8c899caSPeter Brune 
334c8c899caSPeter Brune    Level: advanced
335c8c899caSPeter Brune 
336c8c899caSPeter Brune .keywords: FAS, monitor
337c8c899caSPeter Brune 
338c8c899caSPeter Brune .seealso: SNESFASSetCyclesOnLevel()
339c8c899caSPeter Brune @*/
34022d28d08SBarry Smith PetscErrorCode SNESFASSetMonitor(SNES snes, PetscBool flg)
34122d28d08SBarry Smith {
342c8c899caSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
343c8c899caSPeter Brune   PetscErrorCode ierr;
344c8c899caSPeter Brune   PetscBool      isFine;
345c8c899caSPeter Brune   PetscInt       i, levels = fas->levels;
346c8c899caSPeter Brune   SNES           levelsnes;
34722d28d08SBarry Smith 
348c8c899caSPeter Brune   PetscFunctionBegin;
349c8c899caSPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
350c8c899caSPeter Brune   if (isFine) {
351c8c899caSPeter Brune     for (i = 0; i < levels; i++) {
35222d28d08SBarry Smith       ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
353c8c899caSPeter Brune       fas  = (SNES_FAS*)levelsnes->data;
354c8c899caSPeter Brune       if (flg) {
355ce94432eSBarry Smith         fas->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)levelsnes));CHKERRQ(ierr);
356c8c899caSPeter Brune         /* set the monitors for the upsmoother and downsmoother */
357c8c899caSPeter Brune         ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr);
3580298fd71SBarry Smith         ierr = SNESMonitorSet(levelsnes,SNESMonitorDefault,NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
3591aa26658SKarl Rupp       } else if (i != fas->levels - 1) {
360c8c899caSPeter Brune         /* unset the monitors on the coarse levels */
361c8c899caSPeter Brune         ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr);
362c8c899caSPeter Brune       }
363c8c899caSPeter Brune     }
364c8c899caSPeter Brune   }
365c8c899caSPeter Brune   PetscFunctionReturn(0);
366c8c899caSPeter Brune }
367c8c899caSPeter Brune 
368ab8d36c9SPeter Brune #undef __FUNCT__
3690dd27c6cSPeter Brune #define __FUNCT__ "SNESFASSetLog"
3700dd27c6cSPeter Brune /*@
3710dd27c6cSPeter Brune    SNESFASSetLog - Sets or unsets time logging for various FAS stages on all levels
3720dd27c6cSPeter Brune 
3730dd27c6cSPeter Brune    Logically Collective on SNES
3740dd27c6cSPeter Brune 
3750dd27c6cSPeter Brune    Input Parameters:
3760dd27c6cSPeter Brune +  snes   - the FAS context
3770dd27c6cSPeter Brune -  flg    - monitor or not
3780dd27c6cSPeter Brune 
3790dd27c6cSPeter Brune    Level: advanced
3800dd27c6cSPeter Brune 
3810dd27c6cSPeter Brune .keywords: FAS, logging
3820dd27c6cSPeter Brune 
3830dd27c6cSPeter Brune .seealso: SNESFASSetMonitor()
3840dd27c6cSPeter Brune @*/
3850dd27c6cSPeter Brune PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg)
3860dd27c6cSPeter Brune {
3870dd27c6cSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
3880dd27c6cSPeter Brune   PetscErrorCode ierr;
3890dd27c6cSPeter Brune   PetscBool      isFine;
3900dd27c6cSPeter Brune   PetscInt       i, levels = fas->levels;
3910dd27c6cSPeter Brune   SNES           levelsnes;
3920dd27c6cSPeter Brune   char           eventname[128];
3930dd27c6cSPeter Brune 
3940dd27c6cSPeter Brune   PetscFunctionBegin;
3950dd27c6cSPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
3960dd27c6cSPeter Brune   if (isFine) {
3970dd27c6cSPeter Brune     for (i = 0; i < levels; i++) {
3980dd27c6cSPeter Brune       ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
3990dd27c6cSPeter Brune       fas  = (SNES_FAS*)levelsnes->data;
4000dd27c6cSPeter Brune       if (flg) {
4010dd27c6cSPeter Brune         sprintf(eventname,"FASSetup %d",(int)i);
4020dd27c6cSPeter Brune         ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsetup);CHKERRQ(ierr);
4030dd27c6cSPeter Brune         sprintf(eventname,"FASSmooth %d",(int)i);
4040dd27c6cSPeter Brune         ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsolve);CHKERRQ(ierr);
4050dd27c6cSPeter Brune         sprintf(eventname,"FASResid %d",(int)i);
4060dd27c6cSPeter Brune         ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventresidual);CHKERRQ(ierr);
4070dd27c6cSPeter Brune         if (i) {
4080dd27c6cSPeter Brune           sprintf(eventname,"FASInterp %d",(int)i);
4090dd27c6cSPeter Brune           ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventinterprestrict);CHKERRQ(ierr);
4100dd27c6cSPeter Brune         }
4110dd27c6cSPeter Brune       } else {
4120298fd71SBarry Smith         fas->eventsmoothsetup    = 0;
4130298fd71SBarry Smith         fas->eventsmoothsolve    = 0;
4140298fd71SBarry Smith         fas->eventresidual       = 0;
4150298fd71SBarry Smith         fas->eventinterprestrict = 0;
4160dd27c6cSPeter Brune       }
4170dd27c6cSPeter Brune     }
4180dd27c6cSPeter Brune   }
4190dd27c6cSPeter Brune   PetscFunctionReturn(0);
4200dd27c6cSPeter Brune }
4210dd27c6cSPeter Brune 
4220dd27c6cSPeter Brune #undef __FUNCT__
423ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleCreateSmoother_Private"
424ab8d36c9SPeter Brune /*
425ab8d36c9SPeter Brune Creates the default smoother type.
426ab8d36c9SPeter Brune 
42704d7464bSBarry Smith This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level.
428ab8d36c9SPeter Brune 
429ab8d36c9SPeter Brune  */
43022d28d08SBarry Smith PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth)
43122d28d08SBarry Smith {
432ab8d36c9SPeter Brune   SNES_FAS       *fas;
433ab8d36c9SPeter Brune   const char     *optionsprefix;
434ab8d36c9SPeter Brune   char           tprefix[128];
435ab8d36c9SPeter Brune   PetscErrorCode ierr;
436ab8d36c9SPeter Brune   SNES           nsmooth;
43722d28d08SBarry Smith 
438ab8d36c9SPeter Brune   PetscFunctionBegin;
439ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
440ab8d36c9SPeter Brune   fas  = (SNES_FAS*)snes->data;
441ab8d36c9SPeter Brune   ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr);
442ab8d36c9SPeter Brune   /* create the default smoother */
443ce94432eSBarry Smith   ierr = SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth);CHKERRQ(ierr);
444ab8d36c9SPeter Brune   if (fas->level == 0) {
445ab8d36c9SPeter Brune     sprintf(tprefix,"fas_coarse_");
446ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr);
447ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr);
44804d7464bSBarry Smith     ierr = SNESSetType(nsmooth, SNESNEWTONLS);CHKERRQ(ierr);
449e70c42e5SPeter Brune     ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);CHKERRQ(ierr);
450ab8d36c9SPeter Brune   } else {
451ab8d36c9SPeter Brune     sprintf(tprefix,"fas_levels_%d_",(int)fas->level);
452ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr);
453ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr);
454ab8d36c9SPeter Brune     ierr = SNESSetType(nsmooth, SNESNRICHARDSON);CHKERRQ(ierr);
455e70c42e5SPeter Brune     ierr = SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);CHKERRQ(ierr);
456ab8d36c9SPeter Brune   }
457ab8d36c9SPeter Brune   ierr    = PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
458*3bb1ff40SBarry Smith   ierr    = PetscLogObjectParent((PetscObject)snes,(PetscObject)nsmooth);CHKERRQ(ierr);
459f89ba88eSPeter Brune   ierr    = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth);CHKERRQ(ierr);
460ab8d36c9SPeter Brune   *smooth = nsmooth;
461ab8d36c9SPeter Brune   PetscFunctionReturn(0);
462ab8d36c9SPeter Brune }
463ab8d36c9SPeter Brune 
464ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */
465ab8d36c9SPeter Brune 
466ab8d36c9SPeter Brune #undef __FUNCT__
467ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleSetCycles"
468ab8d36c9SPeter Brune /*@
469ab8d36c9SPeter Brune    SNESFASCycleSetCycles - Sets the number of cycles on a particular level.
470ab8d36c9SPeter Brune 
471ab8d36c9SPeter Brune    Logically Collective on SNES
472ab8d36c9SPeter Brune 
473ab8d36c9SPeter Brune    Input Parameters:
474ab8d36c9SPeter Brune +  snes   - the multigrid context
475ab8d36c9SPeter Brune .  level  - the level to set the number of cycles on
476ab8d36c9SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
477ab8d36c9SPeter Brune 
478ab8d36c9SPeter Brune    Level: advanced
479ab8d36c9SPeter Brune 
480ab8d36c9SPeter Brune .keywords: SNES, FAS, set, cycles, V-cycle, W-cycle, multigrid
481ab8d36c9SPeter Brune 
482ab8d36c9SPeter Brune .seealso: SNESFASSetCycles()
483ab8d36c9SPeter Brune @*/
48422d28d08SBarry Smith PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles)
48522d28d08SBarry Smith {
486ab8d36c9SPeter Brune   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
487ab8d36c9SPeter Brune   PetscErrorCode ierr;
48822d28d08SBarry Smith 
489ab8d36c9SPeter Brune   PetscFunctionBegin;
490ab8d36c9SPeter Brune   fas->n_cycles = cycles;
491ab8d36c9SPeter Brune   ierr          = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr);
492ab8d36c9SPeter Brune   PetscFunctionReturn(0);
493ab8d36c9SPeter Brune }
494ab8d36c9SPeter Brune 
495ab8d36c9SPeter Brune 
496ab8d36c9SPeter Brune #undef __FUNCT__
497ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmoother"
498ab8d36c9SPeter Brune /*@
499ab8d36c9SPeter Brune    SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.
500ab8d36c9SPeter Brune 
501ab8d36c9SPeter Brune    Logically Collective on SNES
502ab8d36c9SPeter Brune 
503ab8d36c9SPeter Brune    Input Parameters:
504ab8d36c9SPeter Brune .  snes   - the multigrid context
505ab8d36c9SPeter Brune 
506ab8d36c9SPeter Brune    Output Parameters:
507ab8d36c9SPeter Brune .  smooth - the smoother
508ab8d36c9SPeter Brune 
509ab8d36c9SPeter Brune    Level: advanced
510ab8d36c9SPeter Brune 
511ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
512ab8d36c9SPeter Brune 
513ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown()
514ab8d36c9SPeter Brune @*/
515ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth)
516ab8d36c9SPeter Brune {
517ab8d36c9SPeter Brune   SNES_FAS *fas;
5185fd66863SKarl Rupp 
519ab8d36c9SPeter Brune   PetscFunctionBegin;
520ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
521ab8d36c9SPeter Brune   fas     = (SNES_FAS*)snes->data;
522ab8d36c9SPeter Brune   *smooth = fas->smoothd;
523ab8d36c9SPeter Brune   PetscFunctionReturn(0);
524ab8d36c9SPeter Brune }
525ab8d36c9SPeter Brune #undef __FUNCT__
526ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmootherUp"
527ab8d36c9SPeter Brune /*@
528ab8d36c9SPeter Brune    SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level.
529ab8d36c9SPeter Brune 
530ab8d36c9SPeter Brune    Logically Collective on SNES
531ab8d36c9SPeter Brune 
532ab8d36c9SPeter Brune    Input Parameters:
533ab8d36c9SPeter Brune .  snes   - the multigrid context
534ab8d36c9SPeter Brune 
535ab8d36c9SPeter Brune    Output Parameters:
536ab8d36c9SPeter Brune .  smoothu - the smoother
537ab8d36c9SPeter Brune 
538ab8d36c9SPeter Brune    Notes:
539ab8d36c9SPeter Brune    Returns the downsmoother if no up smoother is available.  This enables transparent
540ab8d36c9SPeter Brune    default behavior in the process of the solve.
541ab8d36c9SPeter Brune 
542ab8d36c9SPeter Brune    Level: advanced
543ab8d36c9SPeter Brune 
544ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
545ab8d36c9SPeter Brune 
546ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown()
547ab8d36c9SPeter Brune @*/
548ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu)
549ab8d36c9SPeter Brune {
550ab8d36c9SPeter Brune   SNES_FAS *fas;
5515fd66863SKarl Rupp 
552ab8d36c9SPeter Brune   PetscFunctionBegin;
553ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
554ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
5551aa26658SKarl Rupp   if (!fas->smoothu) *smoothu = fas->smoothd;
5561aa26658SKarl Rupp   else *smoothu = fas->smoothu;
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;
5825fd66863SKarl 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:
6050298fd71SBarry Smith    Returns 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;
6165fd66863SKarl 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;
6465fd66863SKarl 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;
6775fd66863SKarl 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;
7085fd66863SKarl 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;
7385fd66863SKarl 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;
7685fd66863SKarl Rupp 
769ab8d36c9SPeter Brune   PetscFunctionBegin;
770ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
771ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
7721aa26658SKarl Rupp   if (fas->level == fas->levels - 1) *flg = PETSC_TRUE;
7731aa26658SKarl Rupp   else *flg = PETSC_FALSE;
774ab8d36c9SPeter Brune   PetscFunctionReturn(0);
775ab8d36c9SPeter Brune }
776ab8d36c9SPeter Brune 
777ab8d36c9SPeter Brune /* ---------- functions called on the finest level that return level-specific information ---------- */
778ab8d36c9SPeter Brune 
779ab8d36c9SPeter Brune #undef __FUNCT__
780ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
781ab8d36c9SPeter Brune /*@
782ab8d36c9SPeter Brune    SNESFASSetInterpolation - Sets the function to be used to calculate the
783ab8d36c9SPeter Brune    interpolation from l-1 to the lth level
784ab8d36c9SPeter Brune 
785ab8d36c9SPeter Brune    Input Parameters:
786ab8d36c9SPeter Brune +  snes      - the multigrid context
787ab8d36c9SPeter Brune .  mat       - the interpolation operator
788ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
789ab8d36c9SPeter Brune 
790ab8d36c9SPeter Brune    Level: advanced
791ab8d36c9SPeter Brune 
792ab8d36c9SPeter Brune    Notes:
793ab8d36c9SPeter Brune           Usually this is the same matrix used also to set the restriction
794ab8d36c9SPeter Brune     for the same level.
795ab8d36c9SPeter Brune 
796ab8d36c9SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
797ab8d36c9SPeter Brune     out from the matrix size which one it is.
798ab8d36c9SPeter Brune 
799ab8d36c9SPeter Brune .keywords:  FAS, multigrid, set, interpolate, level
800ab8d36c9SPeter Brune 
801ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
802ab8d36c9SPeter Brune @*/
8030adebc6cSBarry Smith PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat)
8040adebc6cSBarry Smith {
80522d28d08SBarry Smith   SNES_FAS       *fas;
806ab8d36c9SPeter Brune   PetscErrorCode ierr;
807ab8d36c9SPeter Brune   SNES           levelsnes;
80822d28d08SBarry Smith 
809ab8d36c9SPeter Brune   PetscFunctionBegin;
810ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
811ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
812ab8d36c9SPeter Brune   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
813ab8d36c9SPeter Brune   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
8141aa26658SKarl Rupp 
815ab8d36c9SPeter Brune   fas->interpolate = mat;
816ab8d36c9SPeter Brune   PetscFunctionReturn(0);
817ab8d36c9SPeter Brune }
818ab8d36c9SPeter Brune 
819ab8d36c9SPeter Brune #undef __FUNCT__
820ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetInterpolation"
821ab8d36c9SPeter Brune /*@
822ab8d36c9SPeter Brune    SNESFASGetInterpolation - Gets the matrix used to calculate the
823ab8d36c9SPeter Brune    interpolation from l-1 to the lth level
824ab8d36c9SPeter Brune 
825ab8d36c9SPeter Brune    Input Parameters:
826ab8d36c9SPeter Brune +  snes      - the multigrid context
827ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
828ab8d36c9SPeter Brune 
829ab8d36c9SPeter Brune    Output Parameters:
830ab8d36c9SPeter Brune .  mat       - the interpolation operator
831ab8d36c9SPeter Brune 
832ab8d36c9SPeter Brune    Level: advanced
833ab8d36c9SPeter Brune 
834ab8d36c9SPeter Brune .keywords:  FAS, multigrid, get, interpolate, level
835ab8d36c9SPeter Brune 
836ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale()
837ab8d36c9SPeter Brune @*/
83822d28d08SBarry Smith PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat)
83922d28d08SBarry Smith {
84022d28d08SBarry Smith   SNES_FAS       *fas;
841ab8d36c9SPeter Brune   PetscErrorCode ierr;
842ab8d36c9SPeter Brune   SNES           levelsnes;
84322d28d08SBarry Smith 
844ab8d36c9SPeter Brune   PetscFunctionBegin;
845ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
846ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
847ab8d36c9SPeter Brune   *mat = fas->interpolate;
848ab8d36c9SPeter Brune   PetscFunctionReturn(0);
849ab8d36c9SPeter Brune }
850ab8d36c9SPeter Brune 
851ab8d36c9SPeter Brune #undef __FUNCT__
852ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
853ab8d36c9SPeter Brune /*@
854ab8d36c9SPeter Brune    SNESFASSetRestriction - Sets the function to be used to restrict the defect
855ab8d36c9SPeter Brune    from level l to l-1.
856ab8d36c9SPeter Brune 
857ab8d36c9SPeter Brune    Input Parameters:
858ab8d36c9SPeter Brune +  snes  - the multigrid context
859ab8d36c9SPeter Brune .  mat   - the restriction matrix
860ab8d36c9SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
861ab8d36c9SPeter Brune 
862ab8d36c9SPeter Brune    Level: advanced
863ab8d36c9SPeter Brune 
864ab8d36c9SPeter Brune    Notes:
865ab8d36c9SPeter Brune           Usually this is the same matrix used also to set the interpolation
866ab8d36c9SPeter Brune     for the same level.
867ab8d36c9SPeter Brune 
868ab8d36c9SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
869ab8d36c9SPeter Brune     out from the matrix size which one it is.
870ab8d36c9SPeter Brune 
871ab8d36c9SPeter Brune          If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
872ab8d36c9SPeter Brune     is used.
873ab8d36c9SPeter Brune 
874ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
875ab8d36c9SPeter Brune 
876ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
877ab8d36c9SPeter Brune @*/
87822d28d08SBarry Smith PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat)
87922d28d08SBarry Smith {
88022d28d08SBarry Smith   SNES_FAS       *fas;
881ab8d36c9SPeter Brune   PetscErrorCode ierr;
882ab8d36c9SPeter Brune   SNES           levelsnes;
88322d28d08SBarry Smith 
884ab8d36c9SPeter Brune   PetscFunctionBegin;
885ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
886ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
887ab8d36c9SPeter Brune   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
888ab8d36c9SPeter Brune   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
8891aa26658SKarl Rupp 
890ab8d36c9SPeter Brune   fas->restrct = mat;
891ab8d36c9SPeter Brune   PetscFunctionReturn(0);
892ab8d36c9SPeter Brune }
893ab8d36c9SPeter Brune 
894ab8d36c9SPeter Brune #undef __FUNCT__
895ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetRestriction"
896ab8d36c9SPeter Brune /*@
897ab8d36c9SPeter Brune    SNESFASGetRestriction - Gets the matrix used to calculate the
898ab8d36c9SPeter Brune    restriction from l to the l-1th level
899ab8d36c9SPeter Brune 
900ab8d36c9SPeter Brune    Input Parameters:
901ab8d36c9SPeter Brune +  snes      - the multigrid context
902ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
903ab8d36c9SPeter Brune 
904ab8d36c9SPeter Brune    Output Parameters:
905ab8d36c9SPeter Brune .  mat       - the interpolation operator
906ab8d36c9SPeter Brune 
907ab8d36c9SPeter Brune    Level: advanced
908ab8d36c9SPeter Brune 
909ab8d36c9SPeter Brune .keywords:  FAS, multigrid, get, restrict, level
910ab8d36c9SPeter Brune 
911ab8d36c9SPeter Brune .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale()
912ab8d36c9SPeter Brune @*/
91322d28d08SBarry Smith PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat)
91422d28d08SBarry Smith {
91522d28d08SBarry Smith   SNES_FAS       *fas;
916ab8d36c9SPeter Brune   PetscErrorCode ierr;
917ab8d36c9SPeter Brune   SNES           levelsnes;
91822d28d08SBarry Smith 
919ab8d36c9SPeter Brune   PetscFunctionBegin;
920ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
921ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
922ab8d36c9SPeter Brune   *mat = fas->restrct;
923ab8d36c9SPeter Brune   PetscFunctionReturn(0);
924ab8d36c9SPeter Brune }
925ab8d36c9SPeter Brune 
926ab8d36c9SPeter Brune 
927ab8d36c9SPeter Brune #undef __FUNCT__
928ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInjection"
929ab8d36c9SPeter Brune /*@
930ab8d36c9SPeter Brune    SNESFASSetInjection - Sets the function to be used to inject the solution
931ab8d36c9SPeter Brune    from level l to l-1.
932ab8d36c9SPeter Brune 
933ab8d36c9SPeter Brune    Input Parameters:
934ab8d36c9SPeter Brune  +  snes  - the multigrid context
935ab8d36c9SPeter Brune .  mat   - the restriction matrix
936ab8d36c9SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
937ab8d36c9SPeter Brune 
938ab8d36c9SPeter Brune    Level: advanced
939ab8d36c9SPeter Brune 
940ab8d36c9SPeter Brune    Notes:
941ab8d36c9SPeter Brune          If you do not set this, the restriction and rscale is used to
942ab8d36c9SPeter Brune    project the solution instead.
943ab8d36c9SPeter Brune 
944ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
945ab8d36c9SPeter Brune 
946ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
947ab8d36c9SPeter Brune @*/
94822d28d08SBarry Smith PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat)
94922d28d08SBarry Smith {
95022d28d08SBarry Smith   SNES_FAS       *fas;
951ab8d36c9SPeter Brune   PetscErrorCode ierr;
952ab8d36c9SPeter Brune   SNES           levelsnes;
95322d28d08SBarry Smith 
954ab8d36c9SPeter Brune   PetscFunctionBegin;
955ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
956ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
957ab8d36c9SPeter Brune   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
958ab8d36c9SPeter Brune   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
9591aa26658SKarl Rupp 
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);
10281aa26658SKarl Rupp 
1029ab8d36c9SPeter Brune   fas->rscale = rscale;
1030ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1031ab8d36c9SPeter Brune }
1032ab8d36c9SPeter Brune 
1033ab8d36c9SPeter Brune #undef __FUNCT__
1034ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmoother"
1035ab8d36c9SPeter Brune /*@
1036ab8d36c9SPeter Brune    SNESFASGetSmoother - Gets the default smoother on a level.
1037ab8d36c9SPeter Brune 
1038ab8d36c9SPeter Brune    Input Parameters:
1039ab8d36c9SPeter Brune +  snes   - the multigrid context
1040ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply
1041ab8d36c9SPeter Brune 
1042ab8d36c9SPeter Brune    Output Parameters:
1043ab8d36c9SPeter Brune    smooth  - the smoother
1044ab8d36c9SPeter Brune 
1045ab8d36c9SPeter Brune    Level: advanced
1046ab8d36c9SPeter Brune 
1047ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level
1048ab8d36c9SPeter Brune 
1049ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1050ab8d36c9SPeter Brune @*/
105122d28d08SBarry Smith PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth)
105222d28d08SBarry Smith {
1053ab8d36c9SPeter Brune   SNES_FAS       *fas;
1054ab8d36c9SPeter Brune   PetscErrorCode ierr;
1055ab8d36c9SPeter Brune   SNES           levelsnes;
105622d28d08SBarry Smith 
1057ab8d36c9SPeter Brune   PetscFunctionBegin;
1058ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1059ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
1060ab8d36c9SPeter Brune   if (!fas->smoothd) {
1061ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
1062ab8d36c9SPeter Brune   }
1063ab8d36c9SPeter Brune   *smooth = fas->smoothd;
1064ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1065ab8d36c9SPeter Brune }
1066ab8d36c9SPeter Brune 
1067ab8d36c9SPeter Brune #undef __FUNCT__
1068ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherDown"
1069ab8d36c9SPeter Brune /*@
1070ab8d36c9SPeter Brune    SNESFASGetSmootherDown - Gets the downsmoother on a level.
1071ab8d36c9SPeter Brune 
1072ab8d36c9SPeter Brune    Input Parameters:
1073ab8d36c9SPeter Brune +  snes   - the multigrid context
1074ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply
1075ab8d36c9SPeter Brune 
1076ab8d36c9SPeter Brune    Output Parameters:
1077ab8d36c9SPeter Brune    smooth  - the smoother
1078ab8d36c9SPeter Brune 
1079ab8d36c9SPeter Brune    Level: advanced
1080ab8d36c9SPeter Brune 
1081ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level
1082ab8d36c9SPeter Brune 
1083ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1084ab8d36c9SPeter Brune @*/
108522d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth)
108622d28d08SBarry Smith {
1087ab8d36c9SPeter Brune   SNES_FAS       *fas;
1088ab8d36c9SPeter Brune   PetscErrorCode ierr;
1089ab8d36c9SPeter Brune   SNES           levelsnes;
109022d28d08SBarry Smith 
1091ab8d36c9SPeter Brune   PetscFunctionBegin;
1092ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1093ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
1094ab8d36c9SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1095ab8d36c9SPeter Brune   if (!fas->smoothd) {
1096ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
1097ab8d36c9SPeter Brune   }
1098ab8d36c9SPeter Brune   if (!fas->smoothu) {
1099ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr);
1100ab8d36c9SPeter Brune   }
1101ab8d36c9SPeter Brune   *smooth = fas->smoothd;
1102ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1103ab8d36c9SPeter Brune }
1104ab8d36c9SPeter Brune 
1105ab8d36c9SPeter Brune #undef __FUNCT__
1106ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherUp"
1107ab8d36c9SPeter Brune /*@
1108ab8d36c9SPeter Brune    SNESFASGetSmootherUp - Gets the upsmoother on a level.
1109ab8d36c9SPeter Brune 
1110ab8d36c9SPeter Brune    Input Parameters:
1111ab8d36c9SPeter Brune +  snes   - the multigrid context
1112ab8d36c9SPeter Brune -  level  - the level (0 is coarsest)
1113ab8d36c9SPeter Brune 
1114ab8d36c9SPeter Brune    Output Parameters:
1115ab8d36c9SPeter Brune    smooth  - the smoother
1116ab8d36c9SPeter Brune 
1117ab8d36c9SPeter Brune    Level: advanced
1118ab8d36c9SPeter Brune 
1119ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level
1120ab8d36c9SPeter Brune 
1121ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1122ab8d36c9SPeter Brune @*/
112322d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth)
112422d28d08SBarry Smith {
1125ab8d36c9SPeter Brune   SNES_FAS       *fas;
1126ab8d36c9SPeter Brune   PetscErrorCode ierr;
1127ab8d36c9SPeter Brune   SNES           levelsnes;
112822d28d08SBarry Smith 
1129ab8d36c9SPeter Brune   PetscFunctionBegin;
1130ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1131ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
1132ab8d36c9SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1133ab8d36c9SPeter Brune   if (!fas->smoothd) {
1134ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
1135ab8d36c9SPeter Brune   }
1136ab8d36c9SPeter Brune   if (!fas->smoothu) {
1137ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr);
1138ab8d36c9SPeter Brune   }
1139ab8d36c9SPeter Brune   *smooth = fas->smoothu;
1140ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1141ab8d36c9SPeter Brune }
1142d6ad1212SPeter Brune 
1143d6ad1212SPeter Brune #undef __FUNCT__
1144d6ad1212SPeter Brune #define __FUNCT__ "SNESFASGetCoarseSolve"
1145d6ad1212SPeter Brune /*@
1146d6ad1212SPeter Brune    SNESFASGetCoarseSolve - Gets the coarsest solver.
1147d6ad1212SPeter Brune 
1148d6ad1212SPeter Brune    Input Parameters:
1149d6ad1212SPeter Brune +  snes   - the multigrid context
1150d6ad1212SPeter Brune 
1151d6ad1212SPeter Brune    Output Parameters:
1152d6ad1212SPeter Brune    solve  - the coarse-level solver
1153d6ad1212SPeter Brune 
1154d6ad1212SPeter Brune    Level: advanced
1155d6ad1212SPeter Brune 
1156d6ad1212SPeter Brune .keywords: FAS, MG, get, multigrid, solver, coarse
1157d6ad1212SPeter Brune 
1158d6ad1212SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1159d6ad1212SPeter Brune @*/
116022d28d08SBarry Smith PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *smooth)
116122d28d08SBarry Smith {
1162d6ad1212SPeter Brune   SNES_FAS       *fas;
1163d6ad1212SPeter Brune   PetscErrorCode ierr;
1164d6ad1212SPeter Brune   SNES           levelsnes;
116522d28d08SBarry Smith 
1166d6ad1212SPeter Brune   PetscFunctionBegin;
1167d6ad1212SPeter Brune   ierr = SNESFASGetCycleSNES(snes, 0, &levelsnes);CHKERRQ(ierr);
1168d6ad1212SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
1169d6ad1212SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1170d6ad1212SPeter Brune   if (!fas->smoothd) {
1171d6ad1212SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
1172d6ad1212SPeter Brune   }
1173d6ad1212SPeter Brune   *smooth = fas->smoothd;
1174d6ad1212SPeter Brune   PetscFunctionReturn(0);
1175d6ad1212SPeter Brune }
1176