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