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