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