Level smoothers

From the last tutorial, we have learned that the used multigrid algorithm may have a significant influence in the convergence speed. When comparing the error plots for the standalone multigrid smoothers with unsmoothed and smoothed aggregation multigrid, one finds also a notable difference in the smoothness of the error.

Background on multigrid methods

Obviously, there are cases where some highly oscillatory error modes are left and overlaying some low frequency modes. In other cases there are only low frequency error modes left. Theses are two typical cases one might find in practice.

Multigrid methods are built upont the observation that (cheap) level smoothing method often are able to smooth out high oscillatory error components, whereas they cannot reduce low frequency error components very well. These low frequency error components are then transferred to a coarse level, where they can be seen as high frequency error component for a level smoother on the coarse level.

Warning

Display an image of a sine wave appearing oscillatory on a coarse grid

One should not forget that for an efficient multigrid method both the so-called coarse level correction method and the level smoothers have to work together. That is, one has to choose the right multigrid method (e.g., unsmoothed or sa) in combination with an appropriate level smoothing strategy.

In context of multigrid level smoothers, we have to define both the level smoothers and the coarse solver. Usually, a direct solver is used as coarse solver that is applied to the coarsest multigrid levels. However, it is also possible to apply any other kind of iterative smoothing method or even no solver at all (even though this would be non-standard). The following XML file shows how to use a Jacobi smoother both for level smoothing and as a coarse solver.

Listing 2 ../../../test/tutorial/s1_easy_jacobi.xml
<ParameterList name="MueLu">

  <Parameter name="verbosity" type="string" value="low"/>

  <Parameter name="max levels" type="int" value="10"/>
  <Parameter name="coarse: max size" type="int" value="10"/>

  <Parameter name="multigrid algorithm" type="string" value="unsmoothed"/>

  <!-- Jacobi -->
  <Parameter name="smoother: type" type="string" value="RELAXATION"/>
  <ParameterList name="smoother: params">
    <Parameter name="relaxation: type"  type="string" value="Jacobi"/>
    <Parameter name="relaxation: sweeps" type="int" value="1"/>
    <Parameter name="relaxation: damping factor" type="double" value="0.9"/>
  </ParameterList>

  <!-- Jacobi -->
  <Parameter name="coarse: type" type="string" value="RELAXATION"/>
  <ParameterList name="coarse: params">
    <Parameter name="relaxation: type"  type="string" value="Jacobi"/>
    <Parameter name="relaxation: sweeps" type="int" value="1"/>
    <Parameter name="relaxation: damping factor" type="double" value="0.9"/>
  </ParameterList>

</ParameterList>

The corresponding multigrid hierarchy is

Table 1 and Table 2 show the multigrid effect of different number of Jacobi smoothers on all multigrid levels.

Table 1 2D Laplace equation on 50 x 50 mesh after 1 V-cycle with an AMG multigrid solver and Jacobi smoothers on all multigrid levels (2 processors)
../_images/1level_1jac09.png

Fig. 1 1 level with 1 Jacobi sweep (\(\omega=0.9\))

../_images/1level_10jac09.png

Fig. 2 1 level with 10 Jacobi sweeps (\(\omega=0.9\))

../_images/1level_100jac09.png

Fig. 3 1 level with 100 Jacobi sweeps (\(\omega=0.9\))

../_images/2level_1jac09.png

Fig. 4 2 levels with 1 Jacobi sweep (\(\omega=0.9\))

../_images/2level_10jac09.png

Fig. 5 2 levels with 10 Jacobi sweeps (\(\omega=0.9\))

../_images/2level_100jac09.png

Fig. 6 2 levels with 100 Jacobi sweeps (\(\omega=0.9\))

../_images/3level_1jac09.png

Fig. 7 3 levels with 1 Jacobi sweep (\(\omega=0.9\))

../_images/3level_10jac09.png

Fig. 8 3 levels with 10 Jacobi sweeps (\(\omega=0.9\))

../_images/3level_100jac09.png

Fig. 9 3 levels with 100 Jacobi sweeps (\(\omega=0.9\))

Table 2 2D Laplace equation on 50 x 50 mesh after 5 V-cycle with an AMG multigrid solver and Jacobi smoothers on all multigrid levels. (2 processors)
../_images/5sweeps_1level_1jac09.png

Fig. 10 1 level with 1 Jacobi sweep (\(\omega=0.9\))

../_images/5sweeps_1level_10jac09.png

Fig. 11 1 level with 10 Jacobi sweeps (\(\omega=0.9\))

../_images/5sweeps_1level_100jac09.png

Fig. 12 1 level with 100 Jacobi sweeps (\(\omega=0.9\))

../_images/5sweeps_2level_1jac09.png

Fig. 13 2 levels with 1 Jacobi sweep (\(\omega=0.9\))

../_images/5sweeps_2level_10jac09.png

Fig. 14 2 levels with 10 Jacobi sweeps (\(\omega=0.9\))

../_images/5sweeps_2level_100jac09.png

Fig. 15 2 levels with 100 Jacobi sweeps (\(\omega=0.9\))

../_images/5sweeps_3level_1jac09.png

Fig. 16 3 levels with 1 Jacobi sweep (\(\omega=0.9\))

../_images/5sweeps_3level_10jac09.png

Fig. 17 3 levels with 10 Jacobi sweeps (\(\omega=0.9\))

../_images/5sweeps_3level_100jac09.png

Fig. 18 3 levels with 100 Jacobi sweeps (\(\omega=0.9\))

One has even more fine-grained control over pre- and post-smoothing as shown in the following example, where we use different damping parameters for pre- and post-smoothing (and a direct solver on the coarse grid):

This produces the following multigrid hierarchy

Warning

Insert missing output

Note

Note that the relaxation-based methods provided by the Ifpack/Ifpack2 package are embedded in an outer additive Schwarz method.

Of course, there exist other smoother methods such as polynomial smoothers (Chebyshev) and ILU-based methods. A detailed overview of the different available smoothers can be found in the MueLu User’s Guide ([1]).

Exercise 1

Play around with the smoother parameters and study their effect on the error plot and the convergence of the preconditioned CG method. For all available smoothing options and parameters refer to the MueLu user guide ([1]). Hint: use unsmoothed transfer operator basis functions (i.e., multigrid algorithm = unsmoothed) to highlight the effect of the level smoothers.

Exercise 2

Use the following parameters to solve the \(50\times 50\) Laplace 2D problem on 2 processors:

Listing 3 ../../../test/tutorial/s1_easy_exercise.xml
<ParameterList name="MueLu">

  <Parameter name="verbosity" type="string" value="low"/>

  <Parameter name="max levels" type="int" value="3"/>
  <Parameter name="coarse: max size" type="int" value="10"/>
  <Parameter name="multigrid algorithm" type="string" value="sa"/>

  <!-- Jacobi -->
  <Parameter name="smoother: type" type="string" value="RELAXATION"/>
  <ParameterList name="smoother: params">
    <Parameter name="relaxation: type"  type="string" value="Jacobi"/>
    <Parameter name="relaxation: sweeps" type="int" value="1"/>
    <Parameter name="relaxation: damping factor" type="double" value="0.9"/>
  </ParameterList>

  <!-- Jacobi -->
  <Parameter name="coarse: type" type="string" value="DirectSolver"/>
  
</ParameterList>

That is, we change to smoothed aggregation AMG (SA-AMG). You can find the xml file also in ../../../test/tutorial/s1_easy_exercise.xml. Run the example on two processors and check the number of linear iterations and the solver timings in the screen output. Can you find smoother parameters which reduce the number of iterations? Can you find smoother parameters which reduce the iteration timings?

Footnotes