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