xref: /petsc/src/snes/impls/fas/fas.c (revision a63bb30eff8b1dff6ddce5531d6230e04f6fd10b)
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->next)       ierr = SNESReset(fas->next);CHKERRQ(ierr);
625   PetscFunctionReturn(0);
626 }
627 
628 #undef __FUNCT__
629 #define __FUNCT__ "SNESDestroy_FAS"
630 PetscErrorCode SNESDestroy_FAS(SNES snes)
631 {
632   SNES_FAS * fas = (SNES_FAS *)snes->data;
633   PetscErrorCode ierr = 0;
634 
635   PetscFunctionBegin;
636   /* recursively resets and then destroys */
637   if (fas->upsmooth)     ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
638   if (fas->downsmooth)   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
639   if (fas->inject) {
640     ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
641   }
642   if (fas->interpolate == fas->restrct) {
643     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
644     fas->restrct = PETSC_NULL;
645   } else {
646     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
647     if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
648   }
649   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
650   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
651   ierr = PetscFree(fas);CHKERRQ(ierr);
652 
653   PetscFunctionReturn(0);
654 }
655 
656 #undef __FUNCT__
657 #define __FUNCT__ "SNESSetUp_FAS"
658 PetscErrorCode SNESSetUp_FAS(SNES snes)
659 {
660   SNES_FAS       *fas = (SNES_FAS *) snes->data;
661   PetscErrorCode ierr;
662   VecScatter     injscatter;
663   PetscInt       dm_levels;
664 
665   PetscFunctionBegin;
666 
667   if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) {
668     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
669     dm_levels++;
670     if (dm_levels > fas->levels) {
671       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
672       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
673     }
674   }
675 
676   if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
677     ierr = SNESDefaultGetWork(snes, 1);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
678   } else {
679     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
680   }
681 
682   /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */
683   if (fas->upsmooth) {
684     ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr);
685     if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) {
686       ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
687     }
688     if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) {
689       ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
690     }
691    if (snes->ops->computegs && !fas->upsmooth->ops->computegs) {
692       ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
693     }
694   }
695   if (fas->downsmooth) {
696     ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr);
697     if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) {
698       ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
699     }
700     if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) {
701       ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
702     }
703    if (snes->ops->computegs && !fas->downsmooth->ops->computegs) {
704       ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
705     }
706   }
707   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
708   if (fas->next) {
709     if (fas->galerkin) {
710       ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr);
711     } else {
712       if (snes->ops->computefunction && !fas->next->ops->computefunction) {
713         ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
714       }
715       if (snes->ops->computejacobian && !fas->next->ops->computejacobian) {
716         ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
717       }
718       if (snes->ops->computegs && !fas->next->ops->computegs) {
719         ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
720       }
721     }
722   }
723   if (snes->dm) {
724     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
725     if (fas->next) {
726       /* for now -- assume the DM and the evaluation functions have been set externally */
727       if (!fas->next->dm) {
728         ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
729         ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
730       }
731       /* set the interpolation and restriction from the DM */
732       if (!fas->interpolate) {
733         ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
734         fas->restrct = fas->interpolate;
735       }
736       /* set the injection from the DM */
737       if (!fas->inject) {
738         ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
739         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
740         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
741       }
742     }
743     /* set the DMs of the pre and post-smoothers here */
744     if (fas->upsmooth)  {ierr = SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);}
745     if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);}
746   }
747 
748   /* setup FAS work vectors */
749   if (fas->galerkin) {
750     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
751     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
752   }
753 
754   if (fas->next) {
755    /* gotta set up the solution vector for this to work */
756     if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);}
757     if (!fas->next->vec_rhs) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_rhs);CHKERRQ(ierr);}
758     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
759   }
760   /* got to set them all up at once */
761   PetscFunctionReturn(0);
762 }
763 
764 #undef __FUNCT__
765 #define __FUNCT__ "SNESSetFromOptions_FAS"
766 PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
767 {
768   SNES_FAS       *fas = (SNES_FAS *) snes->data;
769   PetscInt       levels = 1;
770   PetscBool      flg, smoothflg, smoothupflg, smoothdownflg, monflg;
771   PetscErrorCode ierr;
772   const char     *def_smooth = SNESNRICHARDSON;
773   char           pre_type[256];
774   char           post_type[256];
775   char           monfilename[PETSC_MAX_PATH_LEN];
776   SNESFASType    fastype;
777 
778   PetscFunctionBegin;
779   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
780 
781   /* number of levels -- only process on the finest level */
782   if (fas->levels - 1 == fas->level) {
783     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
784     if (!flg && snes->dm) {
785       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
786       levels++;
787       fas->usedmfornumberoflevels = PETSC_TRUE;
788     }
789     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
790   }
791 
792   /* type of pre and/or post smoothers -- set both at once */
793   ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr);
794   ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr);
795   fastype = fas->fastype;
796   ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
797   if (flg) {
798     ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
799   }
800   ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr);
801   if (smoothflg) {
802     ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
803   } else {
804     ierr = PetscOptionsList("-snes_fas_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr);
805     ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr);
806   }
807 
808   /* options for the number of preconditioning cycles and cycle type */
809   ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
810   ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
811   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
812 
813   ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr);
814 
815   ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
816 
817   /* other options for the coarsest level */
818   if (fas->level == 0) {
819     ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr);
820   }
821 
822   ierr = PetscOptionsTail();CHKERRQ(ierr);
823   /* setup from the determined types if there is no pointwise procedure or smoother defined */
824 
825   if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !snes->usegs)) {
826     const char     *prefix;
827     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
828     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
829     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
830     if (fas->level || (fas->levels == 1)) {
831       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_down_");CHKERRQ(ierr);
832     } else {
833       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr);
834     }
835     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
836     ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr);
837   }
838 
839   if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !snes->usegs)) {
840     const char     *prefix;
841     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
842     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
843     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
844     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_up_");CHKERRQ(ierr);
845     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
846     ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr);
847   }
848   if (fas->upsmooth) {
849     ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr);
850   }
851 
852   if (fas->downsmooth) {
853     ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr);
854   }
855 
856   if (fas->level != fas->levels - 1) {
857     ierr = SNESSetTolerances(snes, 0.0, 0.0, 0.0, fas->n_cycles, 1000);CHKERRQ(ierr);
858   }
859 
860   if (monflg) {
861     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
862     /* set the monitors for the upsmoother and downsmoother */
863     ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
864     ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
865     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
866     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
867   } else {
868     /* unset the monitors on the coarse levels */
869     if (fas->level != fas->levels - 1) {
870       ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
871     }
872   }
873 
874   /* recursive option setting for the smoothers */
875   if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);}
876   PetscFunctionReturn(0);
877 }
878 
879 #undef __FUNCT__
880 #define __FUNCT__ "SNESView_FAS"
881 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
882 {
883   SNES_FAS   *fas = (SNES_FAS *) snes->data;
884   PetscBool      iascii;
885   PetscErrorCode ierr;
886 
887   PetscFunctionBegin;
888   if (fas->level != fas->levels - 1) PetscFunctionReturn(0);
889   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
890   if (iascii) {
891     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
892     ierr = PetscViewerASCIIPushTab(viewer);
893     ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
894     if (fas->upsmooth) {
895       ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
896       ierr = PetscViewerASCIIPushTab(viewer);
897       ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
898       ierr = PetscViewerASCIIPopTab(viewer);
899     } else {
900       ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
901     }
902     if (fas->downsmooth) {
903       ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
904       ierr = PetscViewerASCIIPushTab(viewer);
905       ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
906       ierr = PetscViewerASCIIPopTab(viewer);
907     } else {
908       ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
909     }
910     if (snes->usegs) {
911       ierr = PetscViewerASCIIPrintf(viewer, "Using user Gauss-Seidel on level %D -- smoothdown=%D, smoothup=%D\n",
912                                     fas->level, fas->max_down_it, fas->max_up_it);CHKERRQ(ierr);
913     }
914     ierr = PetscViewerASCIIPopTab(viewer);
915   } else {
916     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
917   }
918   PetscFunctionReturn(0);
919 }
920 
921 #undef __FUNCT__
922 #define __FUNCT__ "FASDownSmooth"
923 /*
924 Defines the action of the downsmoother
925  */
926 PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){
927   PetscErrorCode      ierr = 0;
928   PetscReal           fnorm;
929   SNESConvergedReason reason;
930   SNES_FAS            *fas = (SNES_FAS *)snes->data;
931   PetscInt            k;
932   PetscFunctionBegin;
933   if (fas->downsmooth) {
934     ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
935     /* check convergence reason for the smoother */
936     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
937     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
938       snes->reason = SNES_DIVERGED_INNER;
939       PetscFunctionReturn(0);
940     }
941   } else if (snes->usegs && snes->ops->computegs) {
942     if (fas->monitor) {
943       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
944       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
945       ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
946       ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr);
947       ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
948     }
949     for (k = 0; k < fas->max_down_it; k++) {
950       ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
951       if (fas->monitor) {
952         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
953         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
954         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
955         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr);
956         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
957       }
958     }
959   } else if (snes->pc) {
960     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
961     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
962     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
963       snes->reason = SNES_DIVERGED_INNER;
964       PetscFunctionReturn(0);
965     }
966   }
967   PetscFunctionReturn(0);
968 }
969 
970 
971 #undef __FUNCT__
972 #define __FUNCT__ "FASUpSmooth"
973 /*
974 Defines the action of the upsmoother
975  */
976 PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) {
977   PetscErrorCode      ierr = 0;
978   PetscReal           fnorm;
979   SNESConvergedReason reason;
980   SNES_FAS            *fas = (SNES_FAS *)snes->data;
981   PetscInt            k;
982   PetscFunctionBegin;
983   if (fas->upsmooth) {
984     ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
985     /* check convergence reason for the smoother */
986     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
987     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
988       snes->reason = SNES_DIVERGED_INNER;
989       PetscFunctionReturn(0);
990     }
991   } else if (snes->usegs && snes->ops->computegs) {
992     if (fas->monitor) {
993       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
994       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
995       ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
996       ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr);
997       ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
998     }
999     for (k = 0; k < fas->max_up_it; k++) {
1000       ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
1001       if (fas->monitor) {
1002         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1003         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
1004         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
1005         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr);
1006         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
1007       }
1008     }
1009   } else if (snes->pc) {
1010     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
1011     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
1012     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1013       snes->reason = SNES_DIVERGED_INNER;
1014       PetscFunctionReturn(0);
1015     }
1016   }
1017   PetscFunctionReturn(0);
1018 }
1019 
1020 #undef __FUNCT__
1021 #define __FUNCT__ "FASCoarseCorrection"
1022 /*
1023 
1024 Performs the FAS coarse correction as:
1025 
1026 fine problem: F(x) = 0
1027 coarse problem: F^c(x) = b^c
1028 
1029 b^c = F^c(I^c_fx^f - I^c_fF(x))
1030 
1031 with correction:
1032 
1033 
1034 
1035  */
1036 PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
1037   PetscErrorCode      ierr;
1038   Vec                 X_c, Xo_c, F_c, B_c;
1039   SNES_FAS            *fas = (SNES_FAS *)snes->data;
1040   SNESConvergedReason reason;
1041   PetscFunctionBegin;
1042   if (fas->next) {
1043     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1044 
1045     X_c  = fas->next->vec_sol;
1046     Xo_c = fas->next->work[0];
1047     F_c  = fas->next->vec_func;
1048     B_c  = fas->next->vec_rhs;
1049 
1050     /* inject the solution */
1051     if (fas->inject) {
1052       ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr);
1053     } else {
1054       ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
1055       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
1056     }
1057     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
1058 
1059     /* restrict the defect */
1060     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
1061 
1062     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
1063     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
1064     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
1065 
1066     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
1067 
1068     /* set initial guess of the coarse problem to the projected fine solution */
1069     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
1070 
1071     /* recurse to the next level */
1072     fas->next->vec_rhs = B_c;
1073     /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */
1074     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
1075     ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
1076     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1077       snes->reason = SNES_DIVERGED_INNER;
1078       PetscFunctionReturn(0);
1079     }
1080     /* fas->next->vec_rhs = PETSC_NULL; */
1081 
1082     /* correct as x <- x + I(x^c - Rx)*/
1083     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
1084     ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr);
1085   }
1086   PetscFunctionReturn(0);
1087 }
1088 
1089 #undef __FUNCT__
1090 #define __FUNCT__ "FASCycle_Additive"
1091 /*
1092 
1093 The additive cycle looks like:
1094 
1095 xhat = x
1096 xhat = dS(x, b)
1097 x = coarsecorrection(xhat, b_d)
1098 x = x + nu*(xhat - x);
1099 (optional) x = uS(x, b)
1100 
1101 With the coarse RHS (defect correction) as below.
1102 
1103  */
1104 PetscErrorCode FASCycle_Additive(SNES snes, Vec X) {
1105   Vec                 F, B, Xhat;
1106   Vec                 X_c, Xo_c, F_c, B_c, G, W;
1107   PetscErrorCode      ierr;
1108   SNES_FAS *          fas = (SNES_FAS *)snes->data;
1109   SNESConvergedReason reason;
1110   PetscReal           xnorm = 0., fnorm = 0., gnorm = 0., ynorm = 0.;
1111   PetscBool           lssucceed;
1112   PetscFunctionBegin;
1113 
1114   F = snes->vec_func;
1115   B = snes->vec_rhs;
1116   Xhat = snes->work[1];
1117   G    = snes->work[2];
1118   W    = snes->work[3];
1119   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
1120   /* recurse first */
1121   if (fas->next) {
1122     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
1123 
1124     X_c  = fas->next->vec_sol;
1125     Xo_c = fas->next->work[0];
1126     F_c  = fas->next->vec_func;
1127     B_c  = fas->next->vec_rhs;
1128 
1129     /* inject the solution */
1130     if (fas->inject) {
1131       ierr = MatRestrict(fas->inject, Xhat, Xo_c);CHKERRQ(ierr);
1132     } else {
1133       ierr = MatRestrict(fas->restrct, Xhat, Xo_c);CHKERRQ(ierr);
1134       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
1135     }
1136     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
1137 
1138     /* restrict the defect */
1139     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
1140 
1141     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
1142     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
1143     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
1144 
1145     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
1146 
1147     /* set initial guess of the coarse problem to the projected fine solution */
1148     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
1149 
1150     /* recurse */
1151     fas->next->vec_rhs = B_c;
1152     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
1153 
1154     /* smooth on this level */
1155     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
1156 
1157    ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
1158     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1159       snes->reason = SNES_DIVERGED_INNER;
1160       PetscFunctionReturn(0);
1161     }
1162 
1163     /* correct as x <- x + I(x^c - Rx)*/
1164     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
1165     ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr);
1166 
1167     /* additive correction of the coarse direction*/
1168     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1169     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
1170     ierr = VecScale(Xhat, -1.0);CHKERRQ(ierr);
1171     ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Xhat,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
1172     ierr = VecCopy(W, X);CHKERRQ(ierr);
1173     ierr = VecCopy(G, F);CHKERRQ(ierr);
1174     fnorm = gnorm;
1175   } else {
1176     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
1177   }
1178   PetscFunctionReturn(0);
1179 }
1180 
1181 #undef __FUNCT__
1182 #define __FUNCT__ "FASCycle_Multiplicative"
1183 /*
1184 
1185 Defines the FAS cycle as:
1186 
1187 fine problem: F(x) = 0
1188 coarse problem: F^c(x) = b^c
1189 
1190 b^c = F^c(I^c_fx^f - I^c_fF(x))
1191 
1192 correction:
1193 
1194 x = x + I(x^c - Rx)
1195 
1196  */
1197 PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) {
1198 
1199   PetscErrorCode      ierr;
1200   Vec                 F,B;
1201   SNES_FAS            *fas = (SNES_FAS *)snes->data;
1202 
1203   PetscFunctionBegin;
1204   F = snes->vec_func;
1205   B = snes->vec_rhs;
1206   /* pre-smooth -- just update using the pre-smoother */
1207   ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
1208 
1209   ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
1210 
1211   if (fas->level != 0) {
1212     ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr);
1213   }
1214   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1215 
1216   PetscFunctionReturn(0);
1217 }
1218 
1219 #undef __FUNCT__
1220 #define __FUNCT__ "SNESSolve_FAS"
1221 
1222 PetscErrorCode SNESSolve_FAS(SNES snes)
1223 {
1224   PetscErrorCode ierr;
1225   PetscInt       i, maxits;
1226   Vec            X, F;
1227   PetscReal      fnorm;
1228   SNES_FAS       *fas = (SNES_FAS *)snes->data;
1229   PetscFunctionBegin;
1230   maxits = snes->max_its;            /* maximum number of iterations */
1231   snes->reason = SNES_CONVERGED_ITERATING;
1232   X = snes->vec_sol;
1233   F = snes->vec_func;
1234 
1235   /*norm setup */
1236   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1237   snes->iter = 0;
1238   snes->norm = 0.;
1239   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1240   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1241   if (snes->domainerror) {
1242     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1243     PetscFunctionReturn(0);
1244   }
1245   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1246   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1247   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1248   snes->norm = fnorm;
1249   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1250   SNESLogConvHistory(snes,fnorm,0);
1251   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1252 
1253   /* set parameter for default relative tolerance convergence test */
1254   snes->ttol = fnorm*snes->rtol;
1255   /* test convergence */
1256   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1257   if (snes->reason) PetscFunctionReturn(0);
1258   for (i = 0; i < maxits; i++) {
1259     /* Call general purpose update function */
1260 
1261     if (snes->ops->update) {
1262       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
1263     }
1264     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
1265       ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
1266     } else {
1267       ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr);
1268     }
1269 
1270     /* check for FAS cycle divergence */
1271     if (snes->reason != SNES_CONVERGED_ITERATING) {
1272       PetscFunctionReturn(0);
1273     }
1274     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1275     /* Monitor convergence */
1276     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1277     snes->iter = i+1;
1278     snes->norm = fnorm;
1279     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1280     SNESLogConvHistory(snes,snes->norm,0);
1281     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1282     /* Test for convergence */
1283     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1284     if (snes->reason) break;
1285   }
1286   if (i == maxits) {
1287     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
1288     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
1289   }
1290   PetscFunctionReturn(0);
1291 }
1292