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