xref: /petsc/src/snes/impls/fas/fas.c (revision 3dccd265b4a2bee9c2b3da7cb6ddba34177af58c)
1 /* Defines the basic SNES object */
2 #include <../src/snes/impls/fas/fasimpls.h>    /*I  "petscsnesfas.h"  I*/
3 
4 const char *SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","SNESFASType","SNES_FAS",0};
5 
6 /*MC
7 
8 SNESFAS - Full Approximation Scheme nonlinear multigrid solver.
9 
10 The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections
11 
12 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
13 M*/
14 
15 extern PetscErrorCode SNESDestroy_FAS(SNES snes);
16 extern PetscErrorCode SNESSetUp_FAS(SNES snes);
17 extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
18 extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
19 extern PetscErrorCode SNESSolve_FAS(SNES snes);
20 extern PetscErrorCode SNESReset_FAS(SNES snes);
21 extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *);
22 
23 EXTERN_C_BEGIN
24 #undef __FUNCT__
25 #define __FUNCT__ "SNESLineSearchSetType_FAS"
26 PetscErrorCode  SNESLineSearchSetType_FAS(SNES snes, SNESLineSearchType type)
27 {
28   PetscErrorCode ierr;
29   PetscFunctionBegin;
30 
31   switch (type) {
32   case SNES_LS_BASIC:
33     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
34     break;
35   case SNES_LS_BASIC_NONORMS:
36     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
37     break;
38   case SNES_LS_QUADRATIC:
39     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr);
40     break;
41   default:
42     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type.");
43     break;
44   }
45   snes->ls_type = type;
46   PetscFunctionReturn(0);
47 }
48 EXTERN_C_END
49 
50 EXTERN_C_BEGIN
51 
52 #undef __FUNCT__
53 #define __FUNCT__ "SNESCreate_FAS"
54 PetscErrorCode SNESCreate_FAS(SNES snes)
55 {
56   SNES_FAS * fas;
57   PetscErrorCode ierr;
58 
59   PetscFunctionBegin;
60   snes->ops->destroy        = SNESDestroy_FAS;
61   snes->ops->setup          = SNESSetUp_FAS;
62   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
63   snes->ops->view           = SNESView_FAS;
64   snes->ops->solve          = SNESSolve_FAS;
65   snes->ops->reset          = SNESReset_FAS;
66 
67   snes->usesksp             = PETSC_FALSE;
68   snes->usespc              = PETSC_FALSE;
69 
70   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
71   snes->data                  = (void*) fas;
72   fas->level                  = 0;
73   fas->levels                 = 1;
74   fas->n_cycles               = 1;
75   fas->max_up_it              = 1;
76   fas->max_down_it            = 1;
77   fas->upsmooth               = PETSC_NULL;
78   fas->downsmooth             = PETSC_NULL;
79   fas->next                   = PETSC_NULL;
80   fas->previous               = PETSC_NULL;
81   fas->interpolate            = PETSC_NULL;
82   fas->restrct                = PETSC_NULL;
83   fas->inject                 = PETSC_NULL;
84   fas->monitor                = PETSC_NULL;
85   fas->usedmfornumberoflevels = PETSC_FALSE;
86   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
87 
88   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_FAS",SNESLineSearchSetType_FAS);CHKERRQ(ierr);
89   ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr);
90 
91   PetscFunctionReturn(0);
92 }
93 EXTERN_C_END
94 
95 #undef __FUNCT__
96 #define __FUNCT__ "SNESFASGetLevels"
97 /*@
98    SNESFASGetLevels - Gets the number of levels in a FAS.
99 
100    Input Parameter:
101 .  snes - the preconditioner context
102 
103    Output parameter:
104 .  levels - the number of levels
105 
106    Level: advanced
107 
108 .keywords: MG, get, levels, multigrid
109 
110 .seealso: SNESFASSetLevels(), PCMGGetLevels()
111 @*/
112 PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
113   SNES_FAS * fas = (SNES_FAS *)snes->data;
114   PetscFunctionBegin;
115   *levels = fas->levels;
116   PetscFunctionReturn(0);
117 }
118 
119 #undef __FUNCT__
120 #define __FUNCT__ "SNESFASSetCycles"
121 /*@
122    SNESFASSetCycles - Sets the type cycles to use.  Use SNESFASSetCyclesOnLevel() for more
123    complicated cycling.
124 
125    Logically Collective on SNES
126 
127    Input Parameters:
128 +  snes   - the multigrid context
129 -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
130 
131    Options Database Key:
132 $  -snes_fas_cycles 1 or 2
133 
134    Level: advanced
135 
136 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
137 
138 .seealso: SNESFASSetCyclesOnLevel()
139 @*/
140 PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) {
141   SNES_FAS * fas = (SNES_FAS *)snes->data;
142   PetscErrorCode ierr;
143   PetscFunctionBegin;
144   fas->n_cycles = cycles;
145   if (fas->next) {
146     ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr);
147   }
148   PetscFunctionReturn(0);
149 }
150 
151 #undef __FUNCT__
152 #define __FUNCT__ "SNESFASSetCyclesOnLevel"
153 /*@
154    SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level.
155 
156    Logically Collective on SNES
157 
158    Input Parameters:
159 +  snes   - the multigrid context
160 .  level  - the level to set the number of cycles on
161 -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
162 
163    Level: advanced
164 
165 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
166 
167 .seealso: SNESFASSetCycles()
168 @*/
169 PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) {
170   SNES_FAS * fas =  (SNES_FAS *)snes->data;
171   PetscInt top_level = fas->level,i;
172 
173   PetscFunctionBegin;
174   if (level > top_level)
175     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level);
176   /* get to the correct level */
177   for (i = fas->level; i > level; i--) {
178     fas = (SNES_FAS *)fas->next->data;
179   }
180   if (fas->level != level)
181     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel");
182   fas->n_cycles = cycles;
183   PetscFunctionReturn(0);
184 }
185 
186 #undef __FUNCT__
187 #define __FUNCT__ "SNESFASSetGS"
188 /*@C
189    SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used.
190    Use SNESFASSetGSOnLevel() for more complicated staging of smoothers
191    and nonlinear preconditioners.
192 
193    Logically Collective on SNES
194 
195    Input Parameters:
196 +  snes    - the multigrid context
197 .  gsfunc  - the nonlinear smoother function
198 .  ctx     - the user context for the nonlinear smoother
199 -  use_gs  - whether to use the nonlinear smoother or not
200 
201    Level: advanced
202 
203 .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid
204 
205 .seealso: SNESSetGS(), SNESFASSetGSOnLevel()
206 @*/
207 PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
208   PetscErrorCode ierr = 0;
209   SNES_FAS       *fas = (SNES_FAS *)snes->data;
210   PetscFunctionBegin;
211 
212   /* use or don't use it according to user wishes*/
213   snes->usegs = use_gs;
214   if (gsfunc) {
215     ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr);
216     /* push the provided GS up the tree */
217     if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr);
218   } else if (snes->ops->computegs) {
219     /* assume that the user has set the GS solver at this level */
220     if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr);
221   } else if (use_gs) {
222     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level);
223   }
224   PetscFunctionReturn(0);
225 }
226 
227 #undef __FUNCT__
228 #define __FUNCT__ "SNESFASSetGSOnLevel"
229 /*@C
230    SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level.
231 
232    Logically Collective on SNES
233 
234    Input Parameters:
235 +  snes    - the multigrid context
236 .  level   - the level to set the nonlinear smoother on
237 .  gsfunc  - the nonlinear smoother function
238 .  ctx     - the user context for the nonlinear smoother
239 -  use_gs  - whether to use the nonlinear smoother or not
240 
241    Level: advanced
242 
243 .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
244 
245 .seealso: SNESSetGS(), SNESFASSetGS()
246 @*/
247 PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
248   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
249   PetscErrorCode ierr;
250   PetscInt       top_level = fas->level,i;
251   SNES           cur_snes = snes;
252   PetscFunctionBegin;
253   if (level > top_level)
254     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level);
255   /* get to the correct level */
256   for (i = fas->level; i > level; i--) {
257     fas = (SNES_FAS *)fas->next->data;
258     cur_snes = fas->next;
259   }
260   if (fas->level != level)
261     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel");
262   snes->usegs = use_gs;
263   if (gsfunc) {
264     ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr);
265   }
266   PetscFunctionReturn(0);
267 }
268 
269 #undef __FUNCT__
270 #define __FUNCT__ "SNESFASGetSNES"
271 /*@
272    SNESFASGetSNES - Gets the SNES corresponding to a particular
273    level of the FAS hierarchy.
274 
275    Input Parameters:
276 +  snes    - the multigrid context
277    level   - the level to get
278 -  lsnes   - whether to use the nonlinear smoother or not
279 
280    Level: advanced
281 
282 .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
283 
284 .seealso: SNESFASSetLevels(), SNESFASGetLevels()
285 @*/
286 PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
287   SNES_FAS * fas = (SNES_FAS *)snes->data;
288   PetscInt levels = fas->level;
289   PetscInt i;
290   PetscFunctionBegin;
291   *lsnes = snes;
292   if (fas->level < level) {
293     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
294   }
295   if (level > levels - 1) {
296     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
297   }
298   for (i = fas->level; i > level; i--) {
299     *lsnes = fas->next;
300     fas = (SNES_FAS *)(*lsnes)->data;
301   }
302   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
303   PetscFunctionReturn(0);
304 }
305 
306 #undef __FUNCT__
307 #define __FUNCT__ "SNESFASSetType"
308 /*@
309 SNESFASSetType - Sets the update and correction type used for FAS.
310 e
311 
312 
313 @*/
314 PetscErrorCode  SNESFASSetType(SNES snes,SNESFASType fastype)
315 {
316   SNES_FAS                   *fas = (SNES_FAS*)snes->data;
317 
318   PetscFunctionBegin;
319   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
320   PetscValidLogicalCollectiveEnum(snes,fastype,2);
321   fas->fastype = fastype;
322   PetscFunctionReturn(0);
323 }
324 
325 
326 
327 #undef __FUNCT__
328 #define __FUNCT__ "SNESFASSetLevels"
329 /*@C
330    SNESFASSetLevels - Sets the number of levels to use with FAS.
331    Must be called before any other FAS routine.
332 
333    Input Parameters:
334 +  snes   - the snes context
335 .  levels - the number of levels
336 -  comms  - optional communicators for each level; this is to allow solving the coarser
337             problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in
338             Fortran.
339 
340    Level: intermediate
341 
342    Notes:
343      If the number of levels is one then the multigrid uses the -fas_levels prefix
344   for setting the level options rather than the -fas_coarse prefix.
345 
346 .keywords: FAS, MG, set, levels, multigrid
347 
348 .seealso: SNESFASGetLevels()
349 @*/
350 PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
351   PetscErrorCode ierr;
352   PetscInt i;
353   SNES_FAS * fas = (SNES_FAS *)snes->data;
354   SNES prevsnes;
355   MPI_Comm comm;
356   PetscFunctionBegin;
357   comm = ((PetscObject)snes)->comm;
358   if (levels == fas->levels) {
359     if (!comms)
360       PetscFunctionReturn(0);
361   }
362   /* user has changed the number of levels; reset */
363   ierr = SNESReset(snes);CHKERRQ(ierr);
364   /* destroy any coarser levels if necessary */
365   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
366   fas->next = PETSC_NULL;
367   fas->previous = PETSC_NULL;
368   prevsnes = snes;
369   /* setup the finest level */
370   for (i = levels - 1; i >= 0; i--) {
371     if (comms) comm = comms[i];
372     fas->level = i;
373     fas->levels = levels;
374     fas->next = PETSC_NULL;
375     if (i > 0) {
376       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
377       ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr);
378       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
379       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
380       ((SNES_FAS *)fas->next->data)->previous = prevsnes;
381       prevsnes = fas->next;
382       fas = (SNES_FAS *)prevsnes->data;
383     }
384   }
385   PetscFunctionReturn(0);
386 }
387 
388 #undef __FUNCT__
389 #define __FUNCT__ "SNESFASSetNumberSmoothUp"
390 /*@
391    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
392    use on all levels.
393 
394    Logically Collective on PC
395 
396    Input Parameters:
397 +  snes - the multigrid context
398 -  n    - the number of smoothing steps
399 
400    Options Database Key:
401 .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
402 
403    Level: advanced
404 
405 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
406 
407 .seealso: SNESFASSetNumberSmoothDown()
408 @*/
409 PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) {
410   SNES_FAS * fas =  (SNES_FAS *)snes->data;
411   PetscErrorCode ierr = 0;
412   PetscFunctionBegin;
413   fas->max_up_it = n;
414   if (fas->next) {
415     ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr);
416   }
417   PetscFunctionReturn(0);
418 }
419 
420 #undef __FUNCT__
421 #define __FUNCT__ "SNESFASSetNumberSmoothDown"
422 /*@
423    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
424    use on all levels.
425 
426    Logically Collective on PC
427 
428    Input Parameters:
429 +  snes - the multigrid context
430 -  n    - the number of smoothing steps
431 
432    Options Database Key:
433 .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
434 
435    Level: advanced
436 
437 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
438 
439 .seealso: SNESFASSetNumberSmoothUp()
440 @*/
441 PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) {
442   SNES_FAS * fas =  (SNES_FAS *)snes->data;
443   PetscErrorCode ierr = 0;
444   PetscFunctionBegin;
445   fas->max_down_it = n;
446   if (fas->next) {
447     ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr);
448   }
449   PetscFunctionReturn(0);
450 }
451 
452 #undef __FUNCT__
453 #define __FUNCT__ "SNESFASSetInterpolation"
454 /*@
455    SNESFASSetInterpolation - Sets the function to be used to calculate the
456    interpolation from l-1 to the lth level
457 
458    Input Parameters:
459 +  snes      - the multigrid context
460 .  mat       - the interpolation operator
461 -  level     - the level (0 is coarsest) to supply [do not supply 0]
462 
463    Level: advanced
464 
465    Notes:
466           Usually this is the same matrix used also to set the restriction
467     for the same level.
468 
469           One can pass in the interpolation matrix or its transpose; PETSc figures
470     out from the matrix size which one it is.
471 
472 .keywords:  FAS, multigrid, set, interpolate, level
473 
474 .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRscale()
475 @*/
476 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
477   SNES_FAS * fas =  (SNES_FAS *)snes->data;
478   PetscInt top_level = fas->level,i;
479 
480   PetscFunctionBegin;
481   if (level > top_level)
482     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
483   /* get to the correct level */
484   for (i = fas->level; i > level; i--) {
485     fas = (SNES_FAS *)fas->next->data;
486   }
487   if (fas->level != level)
488     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
489   fas->interpolate = mat;
490   PetscFunctionReturn(0);
491 }
492 
493 #undef __FUNCT__
494 #define __FUNCT__ "SNESFASSetRestriction"
495 /*@
496    SNESFASSetRestriction - Sets the function to be used to restrict the defect
497    from level l to l-1.
498 
499    Input Parameters:
500 +  snes  - the multigrid context
501 .  mat   - the restriction matrix
502 -  level - the level (0 is coarsest) to supply [Do not supply 0]
503 
504    Level: advanced
505 
506    Notes:
507           Usually this is the same matrix used also to set the interpolation
508     for the same level.
509 
510           One can pass in the interpolation matrix or its transpose; PETSc figures
511     out from the matrix size which one it is.
512 
513          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
514     is used.
515 
516 .keywords: FAS, MG, set, multigrid, restriction, level
517 
518 .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
519 @*/
520 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
521   SNES_FAS * fas =  (SNES_FAS *)snes->data;
522   PetscInt top_level = fas->level,i;
523 
524   PetscFunctionBegin;
525   if (level > top_level)
526     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
527   /* get to the correct level */
528   for (i = fas->level; i > level; i--) {
529     fas = (SNES_FAS *)fas->next->data;
530   }
531   if (fas->level != level)
532     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
533   fas->restrct = mat;
534   PetscFunctionReturn(0);
535 }
536 
537 
538 #undef __FUNCT__
539 #define __FUNCT__ "SNESFASSetRScale"
540 /*@
541    SNESFASSetRScale - Sets the scaling factor of the restriction
542    operator from level l to l-1.
543 
544    Input Parameters:
545 +  snes   - the multigrid context
546 .  rscale - the restriction scaling
547 -  level  - the level (0 is coarsest) to supply [Do not supply 0]
548 
549    Level: advanced
550 
551    Notes:
552          This is only used in the case that the injection is not set.
553 
554 .keywords: FAS, MG, set, multigrid, restriction, level
555 
556 .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
557 @*/
558 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
559   SNES_FAS * fas =  (SNES_FAS *)snes->data;
560   PetscInt top_level = fas->level,i;
561 
562   PetscFunctionBegin;
563   if (level > top_level)
564     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
565   /* get to the correct level */
566   for (i = fas->level; i > level; i--) {
567     fas = (SNES_FAS *)fas->next->data;
568   }
569   if (fas->level != level)
570     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
571   fas->rscale = rscale;
572   PetscFunctionReturn(0);
573 }
574 
575 
576 #undef __FUNCT__
577 #define __FUNCT__ "SNESFASSetInjection"
578 /*@
579    SNESFASSetInjection - Sets the function to be used to inject the solution
580    from level l to l-1.
581 
582    Input Parameters:
583 +  snes  - the multigrid context
584 .  mat   - the restriction matrix
585 -  level - the level (0 is coarsest) to supply [Do not supply 0]
586 
587    Level: advanced
588 
589    Notes:
590          If you do not set this, the restriction and rscale is used to
591    project the solution instead.
592 
593 .keywords: FAS, MG, set, multigrid, restriction, level
594 
595 .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
596 @*/
597 PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
598   SNES_FAS * fas =  (SNES_FAS *)snes->data;
599   PetscInt top_level = fas->level,i;
600 
601   PetscFunctionBegin;
602   if (level > top_level)
603     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level);
604   /* get to the correct level */
605   for (i = fas->level; i > level; i--) {
606     fas = (SNES_FAS *)fas->next->data;
607   }
608   if (fas->level != level)
609     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection");
610   fas->inject = mat;
611   PetscFunctionReturn(0);
612 }
613 
614 #undef __FUNCT__
615 #define __FUNCT__ "SNESReset_FAS"
616 PetscErrorCode SNESReset_FAS(SNES snes)
617 {
618   PetscErrorCode ierr = 0;
619   SNES_FAS * fas = (SNES_FAS *)snes->data;
620 
621   PetscFunctionBegin;
622   if (fas->upsmooth)     ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr);
623   if (fas->downsmooth)   ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr);
624   if (fas->inject) {
625     ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
626   }
627   if (fas->interpolate == fas->restrct) {
628     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
629     fas->restrct = PETSC_NULL;
630   } else {
631     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
632     if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
633   }
634   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
635 
636   if (fas->upsmooth)   ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr);
637   if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr);
638   if (fas->next)       ierr = SNESReset(fas->next);CHKERRQ(ierr);
639   PetscFunctionReturn(0);
640 }
641 
642 #undef __FUNCT__
643 #define __FUNCT__ "SNESDestroy_FAS"
644 PetscErrorCode SNESDestroy_FAS(SNES snes)
645 {
646   SNES_FAS * fas = (SNES_FAS *)snes->data;
647   PetscErrorCode ierr = 0;
648 
649   PetscFunctionBegin;
650   /* recursively resets and then destroys */
651   if (fas->upsmooth)     ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
652   if (fas->downsmooth)   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
653   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
654   ierr = PetscFree(fas);CHKERRQ(ierr);
655 
656   PetscFunctionReturn(0);
657 }
658 
659 #undef __FUNCT__
660 #define __FUNCT__ "SNESSetUp_FAS"
661 PetscErrorCode SNESSetUp_FAS(SNES snes)
662 {
663   SNES_FAS       *fas = (SNES_FAS *) snes->data;
664   PetscErrorCode ierr;
665   VecScatter     injscatter;
666   PetscInt       dm_levels;
667   Vec            vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
668 
669   PetscFunctionBegin;
670 
671   if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) {
672     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
673     dm_levels++;
674     if (dm_levels > fas->levels) {
675 
676       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESSetUp_FAS*/
677       vec_sol = snes->vec_sol;
678       vec_func = snes->vec_func;
679       vec_sol_update = snes->vec_sol_update;
680       vec_rhs = snes->vec_rhs;
681       snes->vec_sol = PETSC_NULL;
682       snes->vec_func = PETSC_NULL;
683       snes->vec_sol_update = PETSC_NULL;
684       snes->vec_rhs = PETSC_NULL;
685 
686       /* reset the number of levels */
687       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
688       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
689 
690       snes->vec_sol = vec_sol;
691       snes->vec_func = vec_func;
692       snes->vec_rhs = vec_rhs;
693       snes->vec_sol_update = vec_sol_update;
694     }
695   }
696 
697   if (fas->level != fas->levels - 1) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
698 
699   if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
700     ierr = SNESDefaultGetWork(snes, 1);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
701   } else {
702     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
703   }
704 
705   /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */
706   if (fas->upsmooth) {
707     ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr);
708     if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) {
709       ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
710     }
711     if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) {
712       ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
713     }
714    if (snes->ops->computegs && !fas->upsmooth->ops->computegs) {
715       ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
716     }
717   }
718   if (fas->downsmooth) {
719     ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr);
720     if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) {
721       ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
722     }
723     if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) {
724       ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
725     }
726    if (snes->ops->computegs && !fas->downsmooth->ops->computegs) {
727       ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
728     }
729   }
730   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
731   if (fas->next) {
732     if (fas->galerkin) {
733       ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr);
734     } else {
735       if (snes->ops->computefunction && !fas->next->ops->computefunction) {
736         ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
737       }
738       if (snes->ops->computejacobian && !fas->next->ops->computejacobian) {
739         ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
740       }
741       if (snes->ops->computegs && !fas->next->ops->computegs) {
742         ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
743       }
744     }
745   }
746   if (snes->dm) {
747     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
748     if (fas->next) {
749       /* for now -- assume the DM and the evaluation functions have been set externally */
750       if (!fas->next->dm) {
751         ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
752         ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
753       }
754       /* set the interpolation and restriction from the DM */
755       if (!fas->interpolate) {
756         ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
757         fas->restrct = fas->interpolate;
758       }
759       /* set the injection from the DM */
760       if (!fas->inject) {
761         ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
762         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
763         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
764       }
765     }
766     /* set the DMs of the pre and post-smoothers here */
767     if (fas->upsmooth)  {ierr = SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);}
768     if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);}
769   }
770 
771   /* setup FAS work vectors */
772   if (fas->galerkin) {
773     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
774     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
775   }
776 
777   if (fas->next) {
778    /* gotta set up the solution vector for this to work */
779     if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);}
780     if (!fas->next->vec_rhs) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_rhs);CHKERRQ(ierr);}
781     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
782   }
783   /* got to set them all up at once */
784   PetscFunctionReturn(0);
785 }
786 
787 #undef __FUNCT__
788 #define __FUNCT__ "SNESSetFromOptions_FAS"
789 PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
790 {
791   SNES_FAS       *fas = (SNES_FAS *) snes->data;
792   PetscInt       levels = 1;
793   PetscBool      flg, smoothflg, smoothupflg, smoothdownflg, monflg;
794   PetscErrorCode ierr;
795   const char     *def_smooth = SNESNRICHARDSON;
796   char           pre_type[256];
797   char           post_type[256];
798   char           monfilename[PETSC_MAX_PATH_LEN];
799   SNESFASType    fastype;
800 
801   PetscFunctionBegin;
802   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
803 
804   /* number of levels -- only process on the finest level */
805   if (fas->levels - 1 == fas->level) {
806     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
807     if (!flg && snes->dm) {
808       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
809       levels++;
810       fas->usedmfornumberoflevels = PETSC_TRUE;
811     }
812     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
813   }
814 
815   /* type of pre and/or post smoothers -- set both at once */
816   ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr);
817   ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr);
818   fastype = fas->fastype;
819   ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
820   if (flg) {
821     ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
822   }
823   ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr);
824   if (smoothflg) {
825     ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
826   } else {
827     ierr = PetscOptionsList("-snes_fas_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr);
828     ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr);
829   }
830 
831   /* options for the number of preconditioning cycles and cycle type */
832   ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
833   ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
834   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
835 
836   ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr);
837 
838   ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
839 
840   /* other options for the coarsest level */
841   if (fas->level == 0) {
842     ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr);
843   }
844 
845   ierr = PetscOptionsTail();CHKERRQ(ierr);
846   /* setup from the determined types if there is no pointwise procedure or smoother defined */
847 
848   if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !snes->usegs)) {
849     const char     *prefix;
850     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
851     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
852     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
853     if (fas->level || (fas->levels == 1)) {
854       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_down_");CHKERRQ(ierr);
855     } else {
856       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr);
857     }
858     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
859     ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr);
860   }
861 
862   if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !snes->usegs)) {
863     const char     *prefix;
864     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
865     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
866     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
867     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_up_");CHKERRQ(ierr);
868     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
869     ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr);
870   }
871   if (fas->upsmooth) {
872     ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr);
873   }
874 
875   if (fas->downsmooth) {
876     ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr);
877   }
878 
879   if (fas->level != fas->levels - 1) {
880     ierr = SNESSetTolerances(snes, 0.0, 0.0, 0.0, fas->n_cycles, 1000);CHKERRQ(ierr);
881   }
882 
883   if (monflg) {
884     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
885     /* set the monitors for the upsmoother and downsmoother */
886     ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
887     ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
888     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
889     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
890   } else {
891     /* unset the monitors on the coarse levels */
892     if (fas->level != fas->levels - 1) {
893       ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
894     }
895   }
896 
897   /* recursive option setting for the smoothers */
898   if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);}
899   PetscFunctionReturn(0);
900 }
901 
902 #undef __FUNCT__
903 #define __FUNCT__ "SNESView_FAS"
904 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
905 {
906   SNES_FAS   *fas = (SNES_FAS *) snes->data;
907   PetscBool      iascii;
908   PetscErrorCode ierr;
909 
910   PetscFunctionBegin;
911   if (fas->level != fas->levels - 1) PetscFunctionReturn(0);
912   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
913   if (iascii) {
914     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
915     ierr = PetscViewerASCIIPushTab(viewer);
916     ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
917     if (fas->upsmooth) {
918       ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
919       ierr = PetscViewerASCIIPushTab(viewer);
920       ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
921       ierr = PetscViewerASCIIPopTab(viewer);
922     } else {
923       ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
924     }
925     if (fas->downsmooth) {
926       ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
927       ierr = PetscViewerASCIIPushTab(viewer);
928       ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
929       ierr = PetscViewerASCIIPopTab(viewer);
930     } else {
931       ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
932     }
933     if (snes->usegs) {
934       ierr = PetscViewerASCIIPrintf(viewer, "Using user Gauss-Seidel on level %D -- smoothdown=%D, smoothup=%D\n",
935                                     fas->level, fas->max_down_it, fas->max_up_it);CHKERRQ(ierr);
936     }
937     ierr = PetscViewerASCIIPopTab(viewer);
938   } else {
939     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
940   }
941   PetscFunctionReturn(0);
942 }
943 
944 #undef __FUNCT__
945 #define __FUNCT__ "FASDownSmooth"
946 /*
947 Defines the action of the downsmoother
948  */
949 PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){
950   PetscErrorCode      ierr = 0;
951   PetscReal           fnorm;
952   SNESConvergedReason reason;
953   SNES_FAS            *fas = (SNES_FAS *)snes->data;
954   PetscFunctionBegin;
955   if (fas->downsmooth) {
956     ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
957     /* check convergence reason for the smoother */
958     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
959     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
960       snes->reason = SNES_DIVERGED_INNER;
961       PetscFunctionReturn(0);
962     }
963   } else if (snes->usegs && snes->ops->computegs) {
964     if (fas->monitor) {
965       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
966       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
967       ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
968       ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr);
969       ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
970     }
971     ierr = SNESSetGSSweeps(snes, fas->max_down_it);CHKERRQ(ierr);
972     ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
973     if (fas->monitor) {
974       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
975       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
976       ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
977       ierr = PetscViewerASCIIPrintf(fas->monitor, "1 SNES GS Function norm %14.12e\n", fnorm);CHKERRQ(ierr);
978       ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
979     }
980   } else if (snes->pc) {
981     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
982     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
983     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
984       snes->reason = SNES_DIVERGED_INNER;
985       PetscFunctionReturn(0);
986     }
987   }
988   PetscFunctionReturn(0);
989 }
990 
991 
992 #undef __FUNCT__
993 #define __FUNCT__ "FASUpSmooth"
994 /*
995 Defines the action of the upsmoother
996  */
997 PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) {
998   PetscErrorCode      ierr = 0;
999   PetscReal           fnorm;
1000   SNESConvergedReason reason;
1001   SNES_FAS            *fas = (SNES_FAS *)snes->data;
1002   PetscFunctionBegin;
1003   if (fas->upsmooth) {
1004     ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
1005     /* check convergence reason for the smoother */
1006     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
1007     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1008       snes->reason = SNES_DIVERGED_INNER;
1009       PetscFunctionReturn(0);
1010     }
1011   } else if (snes->usegs && snes->ops->computegs) {
1012     if (fas->monitor) {
1013       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1014       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
1015       ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
1016       ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr);
1017       ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
1018     }
1019     ierr = SNESSetGSSweeps(snes, fas->max_up_it);CHKERRQ(ierr);
1020     ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
1021     if (fas->monitor) {
1022       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1023       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
1024       ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
1025       ierr = PetscViewerASCIIPrintf(fas->monitor, "1 SNES GS Function norm %14.12e\n", fnorm);CHKERRQ(ierr);
1026       ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
1027     }
1028   } else if (snes->pc) {
1029     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
1030     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
1031     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1032       snes->reason = SNES_DIVERGED_INNER;
1033       PetscFunctionReturn(0);
1034     }
1035   }
1036   PetscFunctionReturn(0);
1037 }
1038 
1039 #undef __FUNCT__
1040 #define __FUNCT__ "FASCoarseCorrection"
1041 /*
1042 
1043 Performs the FAS coarse correction as:
1044 
1045 fine problem: F(x) = 0
1046 coarse problem: F^c(x) = b^c
1047 
1048 b^c = F^c(I^c_fx^f - I^c_fF(x))
1049 
1050 with correction:
1051 
1052 
1053 
1054  */
1055 PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
1056   PetscErrorCode      ierr;
1057   Vec                 X_c, Xo_c, F_c, B_c;
1058   SNES_FAS            *fas = (SNES_FAS *)snes->data;
1059   SNESConvergedReason reason;
1060   PetscFunctionBegin;
1061   if (fas->next) {
1062     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1063 
1064     X_c  = fas->next->vec_sol;
1065     Xo_c = fas->next->work[0];
1066     F_c  = fas->next->vec_func;
1067     B_c  = fas->next->vec_rhs;
1068 
1069     /* inject the solution */
1070     if (fas->inject) {
1071       ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr);
1072     } else {
1073       ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
1074       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
1075     }
1076     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
1077 
1078     /* restrict the defect */
1079     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
1080 
1081     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
1082     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
1083     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
1084 
1085     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
1086 
1087     /* set initial guess of the coarse problem to the projected fine solution */
1088     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
1089 
1090     /* recurse to the next level */
1091     fas->next->vec_rhs = B_c;
1092     /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */
1093     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
1094     ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
1095     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1096       snes->reason = SNES_DIVERGED_INNER;
1097       PetscFunctionReturn(0);
1098     }
1099     /* fas->next->vec_rhs = PETSC_NULL; */
1100 
1101     /* correct as x <- x + I(x^c - Rx)*/
1102     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
1103     ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr);
1104   }
1105   PetscFunctionReturn(0);
1106 }
1107 
1108 #undef __FUNCT__
1109 #define __FUNCT__ "FASCycle_Additive"
1110 /*
1111 
1112 The additive cycle looks like:
1113 
1114 xhat = x
1115 xhat = dS(x, b)
1116 x = coarsecorrection(xhat, b_d)
1117 x = x + nu*(xhat - x);
1118 (optional) x = uS(x, b)
1119 
1120 With the coarse RHS (defect correction) as below.
1121 
1122  */
1123 PetscErrorCode FASCycle_Additive(SNES snes, Vec X) {
1124   Vec                 F, B, Xhat;
1125   Vec                 X_c, Xo_c, F_c, B_c, G, W;
1126   PetscErrorCode      ierr;
1127   SNES_FAS *          fas = (SNES_FAS *)snes->data;
1128   SNESConvergedReason reason;
1129   PetscReal           xnorm = 0., fnorm = 0., gnorm = 0., ynorm = 0.;
1130   PetscBool           lssucceed;
1131   PetscFunctionBegin;
1132 
1133   F = snes->vec_func;
1134   B = snes->vec_rhs;
1135   Xhat = snes->work[1];
1136   G    = snes->work[2];
1137   W    = snes->work[3];
1138   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
1139   /* recurse first */
1140   if (fas->next) {
1141     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
1142 
1143     X_c  = fas->next->vec_sol;
1144     Xo_c = fas->next->work[0];
1145     F_c  = fas->next->vec_func;
1146     B_c  = fas->next->vec_rhs;
1147 
1148     /* inject the solution */
1149     if (fas->inject) {
1150       ierr = MatRestrict(fas->inject, Xhat, Xo_c);CHKERRQ(ierr);
1151     } else {
1152       ierr = MatRestrict(fas->restrct, Xhat, Xo_c);CHKERRQ(ierr);
1153       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
1154     }
1155     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
1156 
1157     /* restrict the defect */
1158     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
1159 
1160     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
1161     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
1162     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
1163 
1164     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
1165 
1166     /* set initial guess of the coarse problem to the projected fine solution */
1167     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
1168 
1169     /* recurse */
1170     fas->next->vec_rhs = B_c;
1171     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
1172 
1173     /* smooth on this level */
1174     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
1175 
1176    ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
1177     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1178       snes->reason = SNES_DIVERGED_INNER;
1179       PetscFunctionReturn(0);
1180     }
1181 
1182     /* correct as x <- x + I(x^c - Rx)*/
1183     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
1184     ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr);
1185 
1186     /* additive correction of the coarse direction*/
1187     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1188     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
1189     ierr = VecScale(Xhat, -1.0);CHKERRQ(ierr);
1190     ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Xhat,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
1191     ierr = VecCopy(W, X);CHKERRQ(ierr);
1192     ierr = VecCopy(G, F);CHKERRQ(ierr);
1193     fnorm = gnorm;
1194   } else {
1195     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
1196   }
1197   PetscFunctionReturn(0);
1198 }
1199 
1200 #undef __FUNCT__
1201 #define __FUNCT__ "FASCycle_Multiplicative"
1202 /*
1203 
1204 Defines the FAS cycle as:
1205 
1206 fine problem: F(x) = 0
1207 coarse problem: F^c(x) = b^c
1208 
1209 b^c = F^c(I^c_fx^f - I^c_fF(x))
1210 
1211 correction:
1212 
1213 x = x + I(x^c - Rx)
1214 
1215  */
1216 PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) {
1217 
1218   PetscErrorCode      ierr;
1219   Vec                 F,B;
1220   SNES_FAS            *fas = (SNES_FAS *)snes->data;
1221 
1222   PetscFunctionBegin;
1223   F = snes->vec_func;
1224   B = snes->vec_rhs;
1225   /* pre-smooth -- just update using the pre-smoother */
1226   ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
1227 
1228   ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
1229 
1230   if (fas->level != 0) {
1231     ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr);
1232   }
1233   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1234 
1235   PetscFunctionReturn(0);
1236 }
1237 
1238 #undef __FUNCT__
1239 #define __FUNCT__ "SNESSolve_FAS"
1240 
1241 PetscErrorCode SNESSolve_FAS(SNES snes)
1242 {
1243   PetscErrorCode ierr;
1244   PetscInt       i, maxits;
1245   Vec            X, F;
1246   PetscReal      fnorm;
1247   SNES_FAS       *fas = (SNES_FAS *)snes->data;
1248   PetscFunctionBegin;
1249   maxits = snes->max_its;            /* maximum number of iterations */
1250   snes->reason = SNES_CONVERGED_ITERATING;
1251   X = snes->vec_sol;
1252   F = snes->vec_func;
1253 
1254   /*norm setup */
1255   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1256   snes->iter = 0;
1257   snes->norm = 0.;
1258   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1259   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1260   if (snes->domainerror) {
1261     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1262     PetscFunctionReturn(0);
1263   }
1264   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1265   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1266   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1267   snes->norm = fnorm;
1268   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1269   SNESLogConvHistory(snes,fnorm,0);
1270   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1271 
1272   /* set parameter for default relative tolerance convergence test */
1273   snes->ttol = fnorm*snes->rtol;
1274   /* test convergence */
1275   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1276   if (snes->reason) PetscFunctionReturn(0);
1277   for (i = 0; i < maxits; i++) {
1278     /* Call general purpose update function */
1279 
1280     if (snes->ops->update) {
1281       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
1282     }
1283     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
1284       ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
1285     } else {
1286       ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr);
1287     }
1288 
1289     /* check for FAS cycle divergence */
1290     if (snes->reason != SNES_CONVERGED_ITERATING) {
1291       PetscFunctionReturn(0);
1292     }
1293     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1294     /* Monitor convergence */
1295     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1296     snes->iter = i+1;
1297     snes->norm = fnorm;
1298     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1299     SNESLogConvHistory(snes,snes->norm,0);
1300     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1301     /* Test for convergence */
1302     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1303     if (snes->reason) break;
1304   }
1305   if (i == maxits) {
1306     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
1307     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
1308   }
1309   PetscFunctionReturn(0);
1310 }
1311