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