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