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