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