================================ XML interface for advanced users ================================ This tutorial introduces the more advanced (and more flexible) XML interface that can be used for setting up multigrid hierarchies in MueLu. Again we use the 2D Laplace problem on a :math:`50\times 50` mesh as introduced in :ref:`quick_start/example problem`. That is, in the **hands-on.py** script you have to choose option 0 for the problem type. One-level method ================ Before applying a multigrid method as solver, we start with a simple Jacobi iteration as solver and look at the error. By setting the maximum number of multigrid levels to 1 and using a Jacobi smoother as coarse solver we obtain a pseudo multigrid method which corresponds to a simple Jacobi iteration. .. literalinclude:: ../../../test/tutorial/s2_adv_a.xml :language: xml :caption: The advanced XML format is more hierarchical in the structure. Each XML file in the advanced format consists of two major blocks. First there is a set of "Factory" blocks which describe the building blocks within your multigrid methods. In above example there is only one building block specified for the Jacobi method. Each building block needs a (unique) name. In above example the building block has the name **myJacobi**. It is a factory of type **TrilinosSmoother** and describes a damped Jacobi method as declared by the internal parameters. Later we will see examples for other building blocks describing transfer operators or the aggregation strategy. In the "Factories" list, the user has to declare all building blocks of the multigrid method that are used for the setup. The user cannot specify all building blocks involved in the setup process. MueLu will take care of that and use default building blocks for all parts of the setup process where the user makes no explicit statement. This way the user only has to describe what he explicitly needs. It is not sufficient just to declare some building blocks. One also has to register them in the setup process. This is done in the second part of the XML file. The so-called **Hierarchy** block describes the setup process. First there are some basic multigrid parameters that are well-known from the easy XML interface (cf. the previous tutorials). Then, there is an additional list of parameters **All** which encapsulates the information which factory block is responsible to provide certain data. In above example you can see, that the building block **myJacobi** shall be used both for the level smoother and the coarse solver. Since it is only a 1 level problem it would be sufficient to define a coarse solver only. The name of the parameter list **All** can be chosen by the user. It basically describes the user-specified parts of the setup process for all multigrid levels. In this case we just overwrite the internal default factories both for the level smoother and the coarse solver by our Jacobi smoother. .. note:: Be aware that we can have different parameter list sets for different levels, that is, we can use different factories on certain levels. .. admonition:: Exercise 1 Run the **hands-on.py** script and choose the XML file **../../../test/tutorial//s2_adv_a.xml** which contains above XML parameters. Use only 1 processor and visualize the error for an increasing number of multigrid cycles (e.g. 1, 5, 10, 30, 100). What do you observe? .. admonition:: Exercise 2 Note that the relaxation based smoothers are based on a Schwarz method (see Ifpack documentation). Repeat above steps using 2 processors. What do you observe in the error plots? Multigrid method ================ The next step is to introduce a full multigrid algorithm. First one should increase the number of multigrid levels. Second, we switch to a direct solver on the coarsest level. .. admonition:: Exercise 3 Create your own copy of the **../../../test/tutorial/s2_adv_a.xml** parameter file. Adapt it to obtain a 3 level multigrid method. Check how this affects the error plots. .. admonition:: Exercise 4 Change to a direct solver on the coarsest level. You can do this by using **** in the **Hierarchy** block of the xml file. Check the output of the multigrid hierarchy. Level smoothers =============== Next, we give some building blocks for different types of level smoothers that you can use. Note that all these xml blocks can be put into the **Factories** block of the advanced MueLu XML file format. Then you can use them by adding the corresponding link into the **Hierarchy** block using the name of the parameter block. Standard level smoothers ------------------------ Here is a list of the standard level smoothers and how to define them in the XML format. All these standard smoothers are generated by the **TrilinosSmoother** factory class (set in the **factory** parameter in the xml snippets below). - Chebyshev smoother: .. code-block:: xml - Jacobi smoother: .. code-block:: xml - Gauss-Seidel smoother variants: .. code-block:: xml .. code-block:: xml .. code-block:: xml - ILU smoothers: .. code-block:: xml Above listing shows how to create an ILU(0) level smoother (using Ifpack). Please refer to the Ifpack documentation of all parameters on how to choose, e.g. overlapping. .. note:: There is an inconsistency between **Ifpack** (for the **Epetra** stack) and **Ifpack2** (for the **Tpetra** stack) with respect to the **type** parameter in above listing. Not all types of ILU methods available in **Ifpack** are available in **Ifpack2** and vice versa. **Ifpack2** has an implementation of ILUT (**type = ILUT**), but there ILUT parameters may also be slightly different. Please check and adapt the parameters if you switch between the **Epetra** and **Tpetra** stack. Please refer to the **Ifpack** and **Ifpack2** documentation for the details. For Ifpack2 ILUT you can use, e.g., the following settings .. code-block:: xml In Ifpack the parameters might be .. code-block:: xml .. admonition:: Exercise 5 Pick out one level smoother from above and use them for your problem. Note that you may have to adapt the **relaxation: damping factor** for reasonable results. Level smoothers for the Epetra and Tpetra stack ----------------------------------------------- Generally, the **TrilinosSmoother** factory class tries to provide the same set of standard level smoothers for the **Epetra** and **Tpetra** stack as far as possible. It uses **Ifpack** or **Ifpack2** under the hood. **MueLu** partially tries to internally translate the parameters from **Ifpack2** to **Ifpack**, but this is not always possible and error-prone. The available values for the **type** parameter are **RELAXATION**, **CHEBYSHEV**, **ILUT**, **RILUK** and **ILU**. .. note:: In order to define a multigrid hierarchy that is working for both the **Epetra** and **Tpetra** stack it is recommended to use **CHEBYSHEV** or **RELAXATION** based smoothers. The **TrilinosSmoother** factory may provide support more types of smoothers, but these have to be considered as experimental and often are only available for the **Tpetra** branch. Please check the source code for a list of all available smoothers and options. More advanced smoothing strategies will also introduced later in this guide (e.g. :ref:`advanced_concepts/line-smoothing`). Advanced features ================= MueLu allows full control over the behavior of the multigrid levels. Here, we demonstrate the capabilities of MueLu using the level smoothers. Take a look at the following example XML parameter list .. literalinclude:: ../../../test/tutorial/s2_adv_b.xml :language: xml :caption: You can find the parameters in **../../../test/tutorial/s2_adv_b.xml**. We have one building block **BackwardGaussSeidel** representing the level smoother that we want to use in our multigrid hierarchy. As one can see from the **Hierarchy** block we request a 4 level multigrid method. There are two blocks called **Finest** and **Remaining** describing the behavior of the different multigrid levels. Note the **startLevel** parameter in the block **Remaining**. This parameter is missing in the **Finest** block (where it is assumed to be the default value which is zero). That is, in this example we use the backward Gauss-Seidel method as pre-smoother on the finest level (note the keyword **NoFactory** for **PostSmoother**). The parameter **startLevel=1** in the **Remaining** block means that for level 1 and all coarser levels (unless there is another block with **startLevel > 1**) the building blocks from **Remaining** shall be used for the multigrid setup. That is, on the multigrid levels 1 and 2 the backward Gauss-Seidel method is used for post smoothing only. The corresponding multigrid hierarchy has the form .. warning:: \printScreenOutput{s2_adv_b.txt_3.fragment_3.fragment} -> Also replace in Exercise 6. .. admonition:: Exercise 6 Create an XML file in advanced format which produces the following multigrid layout \printScreenOutput{s2_adv_c.txt_3.fragment_3.fragment}. Hint: create a copy of the file **../../../test/tutorial/s2_adv_b.xml** and extend it accordingly. A possible solution can be found in **../../../test/tutorial/s2_adv_c.xml**.