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