xref: /petsc/src/snes/impls/fas/fas.c (revision f4aed992dd926c4e715d02cabab6ea04d956d070)
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   if (gsfunc) {
213     ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr);
214     /* push the provided GS up the tree */
215     if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr);
216   } else if (snes->ops->computegs) {
217     /* assume that the user has set the GS solver at this level */
218     if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr);
219   } else if (use_gs) {
220     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level);
221   }
222   PetscFunctionReturn(0);
223 }
224 
225 #undef __FUNCT__
226 #define __FUNCT__ "SNESFASSetGSOnLevel"
227 /*@C
228    SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level.
229 
230    Logically Collective on SNES
231 
232    Input Parameters:
233 +  snes    - the multigrid context
234 .  level   - the level to set the nonlinear smoother on
235 .  gsfunc  - the nonlinear smoother function
236 .  ctx     - the user context for the nonlinear smoother
237 -  use_gs  - whether to use the nonlinear smoother or not
238 
239    Level: advanced
240 
241 .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
242 
243 .seealso: SNESSetGS(), SNESFASSetGS()
244 @*/
245 PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
246   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
247   PetscErrorCode ierr;
248   PetscInt       top_level = fas->level,i;
249   SNES           cur_snes = snes;
250   PetscFunctionBegin;
251   if (level > top_level)
252     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level);
253   /* get to the correct level */
254   for (i = fas->level; i > level; i--) {
255     fas = (SNES_FAS *)fas->next->data;
256     cur_snes = fas->next;
257   }
258   if (fas->level != level)
259     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel");
260   if (gsfunc) {
261     ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr);
262   }
263   PetscFunctionReturn(0);
264 }
265 
266 #undef __FUNCT__
267 #define __FUNCT__ "SNESFASGetSNES"
268 /*@
269    SNESFASGetSNES - Gets the SNES corresponding to a particular
270    level of the FAS hierarchy.
271 
272    Input Parameters:
273 +  snes    - the multigrid context
274    level   - the level to get
275 -  lsnes   - whether to use the nonlinear smoother or not
276 
277    Level: advanced
278 
279 .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
280 
281 .seealso: SNESFASSetLevels(), SNESFASGetLevels()
282 @*/
283 PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
284   SNES_FAS * fas = (SNES_FAS *)snes->data;
285   PetscInt levels = fas->level;
286   PetscInt i;
287   PetscFunctionBegin;
288   *lsnes = snes;
289   if (fas->level < level) {
290     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
291   }
292   if (level > levels - 1) {
293     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
294   }
295   for (i = fas->level; i > level; i--) {
296     *lsnes = fas->next;
297     fas = (SNES_FAS *)(*lsnes)->data;
298   }
299   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
300   PetscFunctionReturn(0);
301 }
302 
303 #undef __FUNCT__
304 #define __FUNCT__ "SNESFASSetType"
305 /*@
306 SNESFASSetType - Sets the update and correction type used for FAS.
307 e
308 
309 
310 @*/
311 PetscErrorCode  SNESFASSetType(SNES snes,SNESFASType fastype)
312 {
313   SNES_FAS                   *fas = (SNES_FAS*)snes->data;
314 
315   PetscFunctionBegin;
316   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
317   PetscValidLogicalCollectiveEnum(snes,fastype,2);
318   fas->fastype = fastype;
319   PetscFunctionReturn(0);
320 }
321 
322 
323 
324 #undef __FUNCT__
325 #define __FUNCT__ "SNESFASSetLevels"
326 /*@C
327    SNESFASSetLevels - Sets the number of levels to use with FAS.
328    Must be called before any other FAS routine.
329 
330    Input Parameters:
331 +  snes   - the snes context
332 .  levels - the number of levels
333 -  comms  - optional communicators for each level; this is to allow solving the coarser
334             problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in
335             Fortran.
336 
337    Level: intermediate
338 
339    Notes:
340      If the number of levels is one then the multigrid uses the -fas_levels prefix
341   for setting the level options rather than the -fas_coarse prefix.
342 
343 .keywords: FAS, MG, set, levels, multigrid
344 
345 .seealso: SNESFASGetLevels()
346 @*/
347 PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
348   PetscErrorCode ierr;
349   PetscInt i;
350   SNES_FAS * fas = (SNES_FAS *)snes->data;
351   SNES prevsnes;
352   MPI_Comm comm;
353   PetscFunctionBegin;
354   comm = ((PetscObject)snes)->comm;
355   if (levels == fas->levels) {
356     if (!comms)
357       PetscFunctionReturn(0);
358   }
359   /* user has changed the number of levels; reset */
360   ierr = SNESReset(snes);CHKERRQ(ierr);
361   /* destroy any coarser levels if necessary */
362   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
363   fas->next = PETSC_NULL;
364   fas->previous = PETSC_NULL;
365   prevsnes = snes;
366   /* setup the finest level */
367   for (i = levels - 1; i >= 0; i--) {
368     if (comms) comm = comms[i];
369     fas->level = i;
370     fas->levels = levels;
371     fas->next = PETSC_NULL;
372     if (i > 0) {
373       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
374       ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr);
375       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
376       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
377       ((SNES_FAS *)fas->next->data)->previous = prevsnes;
378       prevsnes = fas->next;
379       fas = (SNES_FAS *)prevsnes->data;
380     }
381   }
382   PetscFunctionReturn(0);
383 }
384 
385 #undef __FUNCT__
386 #define __FUNCT__ "SNESFASSetNumberSmoothUp"
387 /*@
388    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
389    use on all levels.
390 
391    Logically Collective on SNES
392 
393    Input Parameters:
394 +  snes - the multigrid context
395 -  n    - the number of smoothing steps
396 
397    Options Database Key:
398 .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
399 
400    Level: advanced
401 
402 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
403 
404 .seealso: SNESFASSetNumberSmoothDown()
405 @*/
406 PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) {
407   SNES_FAS * fas =  (SNES_FAS *)snes->data;
408   PetscErrorCode ierr = 0;
409   PetscFunctionBegin;
410   fas->max_up_it = n;
411   if (fas->next) {
412     ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr);
413   }
414   PetscFunctionReturn(0);
415 }
416 
417 #undef __FUNCT__
418 #define __FUNCT__ "SNESFASSetNumberSmoothDown"
419 /*@
420    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
421    use on all levels.
422 
423    Logically Collective on SNES
424 
425    Input Parameters:
426 +  snes - the multigrid context
427 -  n    - the number of smoothing steps
428 
429    Options Database Key:
430 .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
431 
432    Level: advanced
433 
434 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
435 
436 .seealso: SNESFASSetNumberSmoothUp()
437 @*/
438 PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) {
439   SNES_FAS * fas =  (SNES_FAS *)snes->data;
440   PetscErrorCode ierr = 0;
441   PetscFunctionBegin;
442   fas->max_down_it = n;
443   if (fas->next) {
444     ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr);
445   }
446   PetscFunctionReturn(0);
447 }
448 
449 #undef __FUNCT__
450 #define __FUNCT__ "SNESFASSetInterpolation"
451 /*@
452    SNESFASSetInterpolation - Sets the function to be used to calculate the
453    interpolation from l-1 to the lth level
454 
455    Input Parameters:
456 +  snes      - the multigrid context
457 .  mat       - the interpolation operator
458 -  level     - the level (0 is coarsest) to supply [do not supply 0]
459 
460    Level: advanced
461 
462    Notes:
463           Usually this is the same matrix used also to set the restriction
464     for the same level.
465 
466           One can pass in the interpolation matrix or its transpose; PETSc figures
467     out from the matrix size which one it is.
468 
469 .keywords:  FAS, multigrid, set, interpolate, level
470 
471 .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
472 @*/
473 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
474   SNES_FAS * fas =  (SNES_FAS *)snes->data;
475   PetscInt top_level = fas->level,i;
476   PetscErrorCode ierr;
477 
478   PetscFunctionBegin;
479   if (level > top_level)
480     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
481   /* get to the correct level */
482   for (i = fas->level; i > level; i--) {
483     fas = (SNES_FAS *)fas->next->data;
484   }
485   if (fas->level != level)
486     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
487   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
488   fas->interpolate = mat;
489   PetscFunctionReturn(0);
490 }
491 
492 #undef __FUNCT__
493 #define __FUNCT__ "SNESFASSetRestriction"
494 /*@
495    SNESFASSetRestriction - Sets the function to be used to restrict the defect
496    from level l to l-1.
497 
498    Input Parameters:
499 +  snes  - the multigrid context
500 .  mat   - the restriction matrix
501 -  level - the level (0 is coarsest) to supply [Do not supply 0]
502 
503    Level: advanced
504 
505    Notes:
506           Usually this is the same matrix used also to set the interpolation
507     for the same level.
508 
509           One can pass in the interpolation matrix or its transpose; PETSc figures
510     out from the matrix size which one it is.
511 
512          If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
513     is used.
514 
515 .keywords: FAS, MG, set, multigrid, restriction, level
516 
517 .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
518 @*/
519 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
520   SNES_FAS * fas =  (SNES_FAS *)snes->data;
521   PetscInt top_level = fas->level,i;
522   PetscErrorCode ierr;
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   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
534   fas->restrct = mat;
535   PetscFunctionReturn(0);
536 }
537 
538 
539 #undef __FUNCT__
540 #define __FUNCT__ "SNESFASSetRScale"
541 /*@
542    SNESFASSetRScale - Sets the scaling factor of the restriction
543    operator from level l to l-1.
544 
545    Input Parameters:
546 +  snes   - the multigrid context
547 .  rscale - the restriction scaling
548 -  level  - the level (0 is coarsest) to supply [Do not supply 0]
549 
550    Level: advanced
551 
552    Notes:
553          This is only used in the case that the injection is not set.
554 
555 .keywords: FAS, MG, set, multigrid, restriction, level
556 
557 .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
558 @*/
559 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
560   SNES_FAS * fas =  (SNES_FAS *)snes->data;
561   PetscInt top_level = fas->level,i;
562   PetscErrorCode ierr;
563 
564   PetscFunctionBegin;
565   if (level > top_level)
566     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
567   /* get to the correct level */
568   for (i = fas->level; i > level; i--) {
569     fas = (SNES_FAS *)fas->next->data;
570   }
571   if (fas->level != level)
572     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
573   ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr);
574   fas->rscale = rscale;
575   PetscFunctionReturn(0);
576 }
577 
578 
579 #undef __FUNCT__
580 #define __FUNCT__ "SNESFASSetInjection"
581 /*@
582    SNESFASSetInjection - Sets the function to be used to inject the solution
583    from level l to l-1.
584 
585    Input Parameters:
586 +  snes  - the multigrid context
587 .  mat   - the restriction matrix
588 -  level - the level (0 is coarsest) to supply [Do not supply 0]
589 
590    Level: advanced
591 
592    Notes:
593          If you do not set this, the restriction and rscale is used to
594    project the solution instead.
595 
596 .keywords: FAS, MG, set, multigrid, restriction, level
597 
598 .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
599 @*/
600 PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
601   SNES_FAS * fas =  (SNES_FAS *)snes->data;
602   PetscInt top_level = fas->level,i;
603   PetscErrorCode ierr;
604 
605   PetscFunctionBegin;
606   if (level > top_level)
607     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level);
608   /* get to the correct level */
609   for (i = fas->level; i > level; i--) {
610     fas = (SNES_FAS *)fas->next->data;
611   }
612   if (fas->level != level)
613     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection");
614   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
615   fas->inject = mat;
616   PetscFunctionReturn(0);
617 }
618 
619 #undef __FUNCT__
620 #define __FUNCT__ "SNESReset_FAS"
621 PetscErrorCode SNESReset_FAS(SNES snes)
622 {
623   PetscErrorCode ierr = 0;
624   SNES_FAS * fas = (SNES_FAS *)snes->data;
625 
626   PetscFunctionBegin;
627   ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
628   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
629   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
630   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
631   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
632   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
633 
634   if (fas->upsmooth)   ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr);
635   if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr);
636   if (fas->next)       ierr = SNESReset(fas->next);CHKERRQ(ierr);
637   PetscFunctionReturn(0);
638 }
639 
640 #undef __FUNCT__
641 #define __FUNCT__ "SNESDestroy_FAS"
642 PetscErrorCode SNESDestroy_FAS(SNES snes)
643 {
644   SNES_FAS * fas = (SNES_FAS *)snes->data;
645   PetscErrorCode ierr = 0;
646 
647   PetscFunctionBegin;
648   /* recursively resets and then destroys */
649   ierr = SNESReset(snes);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   Vec            vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
665 
666   PetscFunctionBegin;
667 
668   if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) {
669     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
670     dm_levels++;
671     if (dm_levels > fas->levels) {
672 
673       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESSetUp_FAS*/
674       vec_sol = snes->vec_sol;
675       vec_func = snes->vec_func;
676       vec_sol_update = snes->vec_sol_update;
677       vec_rhs = snes->vec_rhs;
678       snes->vec_sol = PETSC_NULL;
679       snes->vec_func = PETSC_NULL;
680       snes->vec_sol_update = PETSC_NULL;
681       snes->vec_rhs = PETSC_NULL;
682 
683       /* reset the number of levels */
684       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
685       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
686 
687       snes->vec_sol = vec_sol;
688       snes->vec_func = vec_func;
689       snes->vec_rhs = vec_rhs;
690       snes->vec_sol_update = vec_sol_update;
691     }
692   }
693 
694   if (fas->level != fas->levels - 1) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
695 
696   if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
697     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
698   } else {
699     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
700   }
701 
702   if (snes->dm) {
703     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
704     if (fas->next) {
705       /* for now -- assume the DM and the evaluation functions have been set externally */
706       if (!fas->next->dm) {
707         ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
708         ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
709       }
710       /* set the interpolation and restriction from the DM */
711       if (!fas->interpolate) {
712         ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
713         if (!fas->restrct) {
714           ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
715           fas->restrct = fas->interpolate;
716         }
717       }
718       /* set the injection from the DM */
719       if (!fas->inject) {
720         ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
721         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
722         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
723       }
724     }
725     /* set the DMs of the pre and post-smoothers here */
726     if (fas->upsmooth)  {ierr = SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);}
727     if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);}
728   }
729   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
730 
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 
747   if (fas->next) {
748     /* gotta set up the solution vector for this to work */
749     if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);}
750     if (!fas->next->vec_rhs) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_rhs);CHKERRQ(ierr);}
751     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
752   }
753 
754   /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */
755   if (fas->upsmooth) {
756     if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) {
757       ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
758     }
759     if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) {
760       ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
761     }
762    if (snes->ops->computegs && !fas->upsmooth->ops->computegs) {
763       ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
764     }
765    ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr);
766   }
767   if (fas->downsmooth) {
768     if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) {
769       ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
770     }
771     if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) {
772       ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
773     }
774    if (snes->ops->computegs && !fas->downsmooth->ops->computegs) {
775      ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
776     }
777     ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr);
778   }
779 
780   /* setup FAS work vectors */
781   if (fas->galerkin) {
782     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
783     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
784   }
785 
786 
787   /* got to set them all up at once */
788   PetscFunctionReturn(0);
789 }
790 
791 #undef __FUNCT__
792 #define __FUNCT__ "SNESSetFromOptions_FAS"
793 PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
794 {
795   SNES_FAS       *fas = (SNES_FAS *) snes->data;
796   PetscInt       levels = 1;
797   PetscBool      flg, smoothflg, smoothupflg, smoothdownflg, smoothcoarseflg = PETSC_FALSE, monflg;
798   PetscErrorCode ierr;
799   char           monfilename[PETSC_MAX_PATH_LEN];
800   SNESFASType    fastype;
801   const char     *optionsprefix;
802   const char     *prefix;
803 
804   PetscFunctionBegin;
805   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
806 
807   /* number of levels -- only process on the finest level */
808   if (fas->levels - 1 == fas->level) {
809     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
810     if (!flg && snes->dm) {
811       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
812       levels++;
813       fas->usedmfornumberoflevels = PETSC_TRUE;
814     }
815     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
816   }
817   fastype = fas->fastype;
818   ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
819   if (flg) {
820     ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
821   }
822 
823   ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
824 
825   /* smoother setup options */
826   ierr = PetscOptionsHasName(optionsprefix, "-fas_snes_type", &smoothflg);CHKERRQ(ierr);
827   ierr = PetscOptionsHasName(optionsprefix, "-fas_up_snes_type", &smoothupflg);CHKERRQ(ierr);
828   ierr = PetscOptionsHasName(optionsprefix, "-fas_down_snes_type", &smoothdownflg);CHKERRQ(ierr);
829   if (fas->level == 0) {
830     ierr = PetscOptionsHasName(optionsprefix, "-fas_coarse_snes_type", &smoothcoarseflg);CHKERRQ(ierr);
831   }
832   /* options for the number of preconditioning cycles and cycle type */
833 
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 
841   ierr = PetscOptionsTail();CHKERRQ(ierr);
842   /* setup from the determined types if there is no pointwise procedure or smoother defined */
843 
844   if ((!fas->downsmooth) && smoothcoarseflg) {
845     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
846     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
847     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
848     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr);
849     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
850   }
851 
852   if ((!fas->downsmooth) && smoothdownflg) {
853     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
854     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
855     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
856     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_down_");CHKERRQ(ierr);
857     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
858   }
859 
860   if ((!fas->upsmooth) && (fas->level != 0) && smoothupflg) {
861     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
862     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
863     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
864     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_up_");CHKERRQ(ierr);
865     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
866   }
867 
868   if ((!fas->downsmooth) && smoothflg) {
869     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
870     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
871     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
872     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_");CHKERRQ(ierr);
873     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
874   }
875 
876   if ((!fas->upsmooth) && (fas->level != 0) && smoothflg) {
877     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
878     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
879     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
880     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_");CHKERRQ(ierr);
881     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
882   }
883 
884   if (fas->upsmooth) {
885     ierr = SNESSetTolerances(fas->upsmooth, fas->upsmooth->abstol, fas->upsmooth->rtol, fas->upsmooth->xtol, 1, 1000);CHKERRQ(ierr);
886   }
887 
888   if (fas->downsmooth) {
889     ierr = SNESSetTolerances(fas->downsmooth, fas->downsmooth->abstol, fas->downsmooth->rtol, fas->downsmooth->xtol, 1, 1000);CHKERRQ(ierr);
890   }
891 
892   if (fas->level != fas->levels - 1) {
893     ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->xtol, fas->n_cycles, 1000);CHKERRQ(ierr);
894   }
895 
896   /* control the simple Richardson smoother that is default if there's no upsmooth or downsmooth */
897   if (!fas->downsmooth) {
898     ierr = PetscOptionsInt("-fas_down_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
899     ierr = PetscOptionsInt("-fas_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
900     if (fas->level == 0) {
901       ierr = PetscOptionsInt("-fas_coarse_snes_max_it","Coarse smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
902     }
903   }
904 
905   if (!fas->upsmooth) {
906     ierr = PetscOptionsInt("-fas_up_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
907     ierr = PetscOptionsInt("-fas_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
908   }
909 
910   if (monflg) {
911     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
912     /* set the monitors for the upsmoother and downsmoother */
913     ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
914     ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
915     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
916     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
917   } else {
918     /* unset the monitors on the coarse levels */
919     if (fas->level != fas->levels - 1) {
920       ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
921     }
922   }
923 
924   /* recursive option setting for the smoothers */
925   if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);}
926   PetscFunctionReturn(0);
927 }
928 
929 #undef __FUNCT__
930 #define __FUNCT__ "SNESView_FAS"
931 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
932 {
933   SNES_FAS   *fas = (SNES_FAS *) snes->data;
934   PetscBool      iascii;
935   PetscErrorCode ierr;
936 
937   PetscFunctionBegin;
938   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
939   if (iascii) {
940     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
941     ierr = PetscViewerASCIIPushTab(viewer);
942     ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
943     if (fas->upsmooth) {
944       ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
945       ierr = PetscViewerASCIIPushTab(viewer);
946       ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
947       ierr = PetscViewerASCIIPopTab(viewer);
948     } else {
949       ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
950     }
951     if (fas->downsmooth) {
952       ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
953       ierr = PetscViewerASCIIPushTab(viewer);
954       ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
955       ierr = PetscViewerASCIIPopTab(viewer);
956     } else {
957       ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
958     }
959     ierr = PetscViewerASCIIPopTab(viewer);
960   } else {
961     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
962   }
963   PetscFunctionReturn(0);
964 }
965 
966 #undef __FUNCT__
967 #define __FUNCT__ "FASDownSmooth"
968 /*
969 Defines the action of the downsmoother
970  */
971 PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){
972   PetscErrorCode      ierr = 0;
973   PetscReal           fnorm, gnorm, ynorm, xnorm = 0.0;
974   SNESConvergedReason reason;
975   SNES_FAS            *fas = (SNES_FAS *)snes->data;
976   Vec                 G, W, Y;
977   PetscBool           lssuccess;
978   PetscInt            k;
979   PetscFunctionBegin;
980   G = snes->work[1];
981   W = snes->work[2];
982   Y = snes->work[3];
983   if (fas->downsmooth) {
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 {
992     /* basic richardson smoother */
993     for (k = 0; k < fas->max_up_it; k++) {
994       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
995       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
996       ierr = VecCopy(F, Y);CHKERRQ(ierr);
997       ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr);
998       if (!lssuccess) {
999         if (++snes->numFailures >= snes->maxFailures) {
1000           snes->reason = SNES_DIVERGED_LINE_SEARCH;
1001           PetscFunctionReturn(0);
1002         }
1003       }
1004       ierr = VecCopy(W, X);CHKERRQ(ierr);
1005       ierr = VecCopy(G, F);CHKERRQ(ierr);
1006       fnorm = gnorm;
1007     }
1008   }
1009   PetscFunctionReturn(0);
1010 }
1011 
1012 
1013 #undef __FUNCT__
1014 #define __FUNCT__ "FASUpSmooth"
1015 /*
1016 Defines the action of the upsmoother
1017  */
1018 PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) {
1019   PetscErrorCode      ierr = 0;
1020   PetscReal           fnorm, gnorm, ynorm, xnorm = 0.0;
1021   SNESConvergedReason reason;
1022   SNES_FAS            *fas = (SNES_FAS *)snes->data;
1023   Vec                 G, W, Y;
1024   PetscBool           lssuccess;
1025   PetscInt            k;
1026   PetscFunctionBegin;
1027   G = snes->work[1];
1028   W = snes->work[2];
1029   Y = snes->work[3];
1030   if (fas->upsmooth) {
1031     ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
1032     /* check convergence reason for the smoother */
1033     ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr);
1034     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1035       snes->reason = SNES_DIVERGED_INNER;
1036       PetscFunctionReturn(0);
1037     }
1038   } else {
1039     /* basic richardson smoother */
1040     for (k = 0; k < fas->max_up_it; k++) {
1041       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1042       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
1043       ierr = VecCopy(F, Y);CHKERRQ(ierr);
1044       ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr);
1045       if (!lssuccess) {
1046         if (++snes->numFailures >= snes->maxFailures) {
1047           snes->reason = SNES_DIVERGED_LINE_SEARCH;
1048           PetscFunctionReturn(0);
1049         }
1050       }
1051       ierr = VecCopy(W, X);CHKERRQ(ierr);
1052       ierr = VecCopy(G, F);CHKERRQ(ierr);
1053       fnorm = gnorm;
1054     }
1055   }
1056   PetscFunctionReturn(0);
1057 }
1058 
1059 #undef __FUNCT__
1060 #define __FUNCT__ "FASCoarseCorrection"
1061 /*
1062 
1063 Performs the FAS coarse correction as:
1064 
1065 fine problem: F(x) = 0
1066 coarse problem: F^c(x) = b^c
1067 
1068 b^c = F^c(I^c_fx^f - I^c_fF(x))
1069 
1070 with correction:
1071 
1072 
1073 
1074  */
1075 PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
1076   PetscErrorCode      ierr;
1077   Vec                 X_c, Xo_c, F_c, B_c;
1078   SNES_FAS            *fas = (SNES_FAS *)snes->data;
1079   SNESConvergedReason reason;
1080   PetscFunctionBegin;
1081   if (fas->next) {
1082     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1083 
1084     X_c  = fas->next->vec_sol;
1085     Xo_c = fas->next->work[0];
1086     F_c  = fas->next->vec_func;
1087     B_c  = fas->next->vec_rhs;
1088 
1089     /* inject the solution */
1090     if (fas->inject) {
1091       ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr);
1092     } else {
1093       ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
1094       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
1095     }
1096     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
1097 
1098     /* restrict the defect */
1099     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
1100 
1101     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
1102     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
1103     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
1104 
1105     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
1106 
1107     /* set initial guess of the coarse problem to the projected fine solution */
1108     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
1109 
1110     /* recurse to the next level */
1111     fas->next->vec_rhs = B_c;
1112     /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */
1113     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
1114     ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
1115     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1116       snes->reason = SNES_DIVERGED_INNER;
1117       PetscFunctionReturn(0);
1118     }
1119     /* fas->next->vec_rhs = PETSC_NULL; */
1120 
1121     /* correct as x <- x + I(x^c - Rx)*/
1122     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
1123     ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr);
1124   }
1125   PetscFunctionReturn(0);
1126 }
1127 
1128 #undef __FUNCT__
1129 #define __FUNCT__ "FASCycle_Additive"
1130 /*
1131 
1132 The additive cycle looks like:
1133 
1134 xhat = x
1135 xhat = dS(x, b)
1136 x = coarsecorrection(xhat, b_d)
1137 x = x + nu*(xhat - x);
1138 (optional) x = uS(x, b)
1139 
1140 With the coarse RHS (defect correction) as below.
1141 
1142  */
1143 PetscErrorCode FASCycle_Additive(SNES snes, Vec X) {
1144   Vec                 F, B, Xhat;
1145   Vec                 X_c, Xo_c, F_c, B_c, G, W;
1146   PetscErrorCode      ierr;
1147   SNES_FAS *          fas = (SNES_FAS *)snes->data;
1148   SNESConvergedReason reason;
1149   PetscReal           xnorm = 0., fnorm = 0., gnorm = 0., ynorm = 0.;
1150   PetscBool           lssucceed;
1151   PetscFunctionBegin;
1152 
1153   F = snes->vec_func;
1154   B = snes->vec_rhs;
1155   Xhat = snes->work[3];
1156   G    = snes->work[1];
1157   W    = snes->work[2];
1158   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
1159   /* recurse first */
1160   if (fas->next) {
1161     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
1162 
1163     X_c  = fas->next->vec_sol;
1164     Xo_c = fas->next->work[0];
1165     F_c  = fas->next->vec_func;
1166     B_c  = fas->next->vec_rhs;
1167 
1168     /* inject the solution */
1169     if (fas->inject) {
1170       ierr = MatRestrict(fas->inject, Xhat, Xo_c);CHKERRQ(ierr);
1171     } else {
1172       ierr = MatRestrict(fas->restrct, Xhat, Xo_c);CHKERRQ(ierr);
1173       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
1174     }
1175     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
1176 
1177     /* restrict the defect */
1178     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
1179 
1180     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
1181     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
1182     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
1183 
1184     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
1185 
1186     /* set initial guess of the coarse problem to the projected fine solution */
1187     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
1188 
1189     /* recurse */
1190     fas->next->vec_rhs = B_c;
1191     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
1192 
1193     /* smooth on this level */
1194     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
1195 
1196    ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
1197     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1198       snes->reason = SNES_DIVERGED_INNER;
1199       PetscFunctionReturn(0);
1200     }
1201 
1202     /* correct as x <- x + I(x^c - Rx)*/
1203     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
1204     ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr);
1205 
1206     /* additive correction of the coarse direction*/
1207     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1208     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
1209     ierr = VecScale(Xhat, -1.0);CHKERRQ(ierr);
1210     ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Xhat,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
1211     ierr = VecCopy(W, X);CHKERRQ(ierr);
1212     ierr = VecCopy(G, F);CHKERRQ(ierr);
1213     fnorm = gnorm;
1214   } else {
1215     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
1216   }
1217   PetscFunctionReturn(0);
1218 }
1219 
1220 #undef __FUNCT__
1221 #define __FUNCT__ "FASCycle_Multiplicative"
1222 /*
1223 
1224 Defines the FAS cycle as:
1225 
1226 fine problem: F(x) = 0
1227 coarse problem: F^c(x) = b^c
1228 
1229 b^c = F^c(I^c_fx^f - I^c_fF(x))
1230 
1231 correction:
1232 
1233 x = x + I(x^c - Rx)
1234 
1235  */
1236 PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) {
1237 
1238   PetscErrorCode      ierr;
1239   Vec                 F,B;
1240   SNES_FAS            *fas = (SNES_FAS *)snes->data;
1241 
1242   PetscFunctionBegin;
1243   F = snes->vec_func;
1244   B = snes->vec_rhs;
1245   /* pre-smooth -- just update using the pre-smoother */
1246   ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
1247 
1248   ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
1249 
1250   if (fas->level != 0) {
1251     ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr);
1252   }
1253   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1254 
1255   PetscFunctionReturn(0);
1256 }
1257 
1258 #undef __FUNCT__
1259 #define __FUNCT__ "SNESSolve_FAS"
1260 
1261 PetscErrorCode SNESSolve_FAS(SNES snes)
1262 {
1263   PetscErrorCode ierr;
1264   PetscInt       i, maxits;
1265   Vec            X, F;
1266   PetscReal      fnorm;
1267   SNES_FAS       *fas = (SNES_FAS *)snes->data;
1268   PetscFunctionBegin;
1269   maxits = snes->max_its;            /* maximum number of iterations */
1270   snes->reason = SNES_CONVERGED_ITERATING;
1271   X = snes->vec_sol;
1272   F = snes->vec_func;
1273 
1274   /*norm setup */
1275   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1276   snes->iter = 0;
1277   snes->norm = 0.;
1278   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1279   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1280   if (snes->domainerror) {
1281     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1282     PetscFunctionReturn(0);
1283   }
1284   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1285   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1286   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1287   snes->norm = fnorm;
1288   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1289   SNESLogConvHistory(snes,fnorm,0);
1290   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1291 
1292   /* set parameter for default relative tolerance convergence test */
1293   snes->ttol = fnorm*snes->rtol;
1294   /* test convergence */
1295   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1296   if (snes->reason) PetscFunctionReturn(0);
1297   for (i = 0; i < maxits; i++) {
1298     /* Call general purpose update function */
1299 
1300     if (snes->ops->update) {
1301       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
1302     }
1303     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
1304       ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
1305     } else {
1306       ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr);
1307     }
1308 
1309     /* check for FAS cycle divergence */
1310     if (snes->reason != SNES_CONVERGED_ITERATING) {
1311       PetscFunctionReturn(0);
1312     }
1313     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1314     /* Monitor convergence */
1315     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1316     snes->iter = i+1;
1317     snes->norm = fnorm;
1318     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1319     SNESLogConvHistory(snes,snes->norm,0);
1320     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1321     /* Test for convergence */
1322     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1323     if (snes->reason) break;
1324   }
1325   if (i == maxits) {
1326     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
1327     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
1328   }
1329   PetscFunctionReturn(0);
1330 }
1331