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