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