Rebalancing

Basic concepts and parameters

Especially when using the uncoupled aggregation strategy it is essential to reduce the number of processors used on the coarser levels. A natural strategy is to make sure that each processor gets a minimum number of equations to solve and reduce the number of active processors accordingly.

In this tutorial we use a multijagged-based repartitioning for the coarse level matrices $A_c$ to rebalance the coarse level problems. The repartitioning algorithm is implemented in the Zoltan2 package of Trilinos. A potential disadvantage of the multijagged-based repartitioning methods is that the algorithm needs additional node coordinates. There are algorithms such as hypergraph partitioning that don’t require coordinates, but current implementations are much more computationally expensive without a commensurate increase in quality.

Repartitioning algorithms are a very wide field of research and can be very complicated. Here, we cannot go into details and just focus on how to use them. Basically there are only a few really important parameters that the user has to set properly:

Description

repartition: min rows per proc

The minimum number of rows each processor shall handle. This parameter is used to reduce the number of involved processors on the coarser levels. If for example the parameter value is chosen to be 1000 and the fine level problem has 10000 rows whereas the coarse level problem has 2000 rows, then the fine level problem is solved on not more than 10 processors at maximum and for the coarse level problem there are not more processors than at maximum 2 being used.

repartition: max imbalance

This parameter defines the maximum allowed imbalance ratio of nonzeros on all processors. If the value is set to 1.2, and there is one processor with more than 20% nonzeros compared to another processor, then the problem will be rebalanced.

repartition: start level

Start rebalancing on the given level and coarser levels. This allows to avoid the costs of rebalancing on the finer levels (where it is not really necessary).

Transfer operator design

Warning

Insert and link missing figures

Figure ref{fig:rebalanceddesignpgamg} gives the extended factory design for smoothed aggregation based AMG for non-symmetric linear systems with rebalancing enabled. Nothing has changed in the upper part where the non-rebalanced Galerkin product has been calculated using the RAPFactory. The coarse level matrix \(A_c\) as output from the RAPFactory then is checked for its partition and rebalanced.

The AmalgamationFactory amalgamates the matrix, i.e. it generates some mapping between the actual degrees of freedom and the corresponding nodes or supernodes. In fact the AmalgamationFactory is only important if there are more than one degree of freedom per node. Otherwise the mappings are trivial to build.

The RepartitionHeuristicFactory contains the rebalancing logic. Depending on the chosen repartitioning parameters, it determines the number of partitions (variable name number of partitions) for the coarse level problem. The IsorropiaInterface class first builds internally the graph of the coarse level matrix \(A_c\) using the information from the AmalgamationInformation and then calls the repartitioning algorithm from Zoltan through the Isorropia interface. The output is an amalgamated repartitioning information. Then the RepartitionInterface factory resembles the un-amalgamated repartitioning information which is put into the RepartitionFactory.

The RepartitionFactory creates the communication “plan” that is used to rebalance the transfer operators and the coarse level matrix.

Insert missing figure here

XML interface

The corresponding XML parameter file looks like

Listing 7 ../../../test/tutorial/s5a.xml
<ParameterList name="MueLu">

  <!-- Factory collection -->
  <ParameterList name="Factories">

    <ParameterList name="myTentativePFact">
      <Parameter name="factory"                             type="string" value="TentativePFactory"/>
    </ParameterList>
    <ParameterList name="myProlongatorFact">
      <Parameter name="factory"                             type="string" value="PgPFactory"/>
      <Parameter name="P"                                   type="string" value="myTentativePFact"/>
    </ParameterList>
    <ParameterList name="myRestrictorFact">
      <Parameter name="factory"                             type="string" value="GenericRFactory"/>
      <Parameter name="P"                                   type="string" value="myProlongatorFact"/>
    </ParameterList>

    <ParameterList name="myAggExportFact">
      <Parameter name="factory"                             type="string" value="AggregationExportFactory"/>
      <Parameter name="Output filename"                     type="string" value="aggs_level%LEVELID_proc%PROCID.out"/>
    </ParameterList>

    <ParameterList name="myRAPFact">
      <Parameter name="factory"                             type="string" value="RAPFactory"/>
      <Parameter name="P"                                   type="string" value="myProlongatorFact"/>
      <Parameter name="R"                                   type="string" value="myRestrictorFact"/>
      <ParameterList name="TransferFactories">
        <Parameter name="Visualization"                     type="string" value="myAggExportFact"/>
      </ParameterList>
    </ParameterList>

    <!-- =======================  REPARTITIONING  ======================= -->
    <!-- amalgamation of coarse level matrix -->
    <ParameterList name="myRebAmalgFact">
      <Parameter name="factory"                        type="string" value="AmalgamationFactory"/>
      <Parameter name="A"                              type="string" value="myRAPFact"/>
    </ParameterList>

    <ParameterList name="myRepartitionHeuristicFact">
      <Parameter name="factory"                        type="string" value="RepartitionHeuristicFactory"/>
      <Parameter name="A"                              type="string" value="myRAPFact"/>
      <Parameter name="repartition: min rows per proc"      type="int"    value="2000"/>
      <Parameter name="repartition: max imbalance"          type="double" value="1.1"/>
      <Parameter name="repartition: start level"            type="int"    value="1"/>    
    </ParameterList>
    
    <ParameterList name="myIsorropiaInterface">
      <Parameter name="factory"                        type="string" value="IsorropiaInterface"/>
      <Parameter name="A"                              type="string" value="myRAPFact"/>
      <Parameter name="UnAmalgamationInfo"             type="string" value="myRebAmalgFact"/>
      <Parameter name="number of partitions"           type="string" value="myRepartitionHeuristicFact"/>
    </ParameterList>

    <ParameterList name="myRepartitionInterface">
      <Parameter name="factory"                        type="string" value="RepartitionInterface"/>
      <Parameter name="A"                              type="string" value="myRAPFact"/>
      <Parameter name="AmalgamatedPartition"           type="string" value="myIsorropiaInterface"/>
      <Parameter name="number of partitions"           type="string" value="myRepartitionHeuristicFact"/>
    </ParameterList>


    <ParameterList name="myRepartitionFact">
      <Parameter name="factory"                        type="string" value="RepartitionFactory"/>
      <Parameter name="A"                              type="string" value="myRAPFact"/>
      <Parameter name="Partition"                      type="string" value="myRepartitionInterface"/>
      <Parameter name="repartition: remap parts"       type="bool"   value="false"/>
      <Parameter name="number of partitions"           type="string" value="myRepartitionHeuristicFact"/>
    </ParameterList>

    <ParameterList name="myRebalanceProlongatorFact">
      <Parameter name="factory"                        type="string" value="RebalanceTransferFactory"/>
      <Parameter name="type"                           type="string" value="Interpolation"/>
      <Parameter name="P"                              type="string" value="myProlongatorFact"/>
      <Parameter name="Nullspace"                      type="string" value="myTentativePFact"/>
    </ParameterList>

    <ParameterList name="myRebalanceRestrictionFact">
      <Parameter name="factory"                        type="string" value="RebalanceTransferFactory"/>
      <Parameter name="type"                           type="string" value="Restriction"/>
      <Parameter name="R"                              type="string" value="myRestrictorFact"/>
    </ParameterList>

    <ParameterList name="myRebalanceAFact">
      <Parameter name="factory"                        type="string" value="RebalanceAcFactory"/>
      <Parameter name="A"                              type="string" value="myRAPFact"/>
    </ParameterList>

    <!-- =======================  SMOOTHERS  ======================= -->
    <ParameterList name="SymGaussSeidel">
      <Parameter name="factory"                        type="string" value="TrilinosSmoother"/>
      <Parameter name="type"                           type="string" value="RELAXATION"/>
      <ParameterList name="ParameterList">
        <Parameter name="relaxation: type"             type="string" value="Symmetric Gauss-Seidel"/>
        <Parameter name="relaxation: sweeps"           type="int"    value="20"/>
        <Parameter name="relaxation: damping factor"   type="double" value="1.2"/>
      </ParameterList>
    </ParameterList>

  </ParameterList>

  <!-- Definition of the multigrid preconditioner -->
  <ParameterList name="Hierarchy">

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

    <ParameterList name="All">
      <Parameter name="Smoother"                            type="string"   value="SymGaussSeidel"/>
      <Parameter name="Nullspace"                           type="string"   value="myRebalanceProlongatorFact"/>
      <Parameter name="P"                                   type="string"   value="myRebalanceProlongatorFact"/>
      <Parameter name="R"                                   type="string"   value="myRebalanceRestrictionFact"/>
      <Parameter name="A"                                   type="string"   value="myRebalanceAFact"/>
      <Parameter name="Importer"                            type="string"   value="myRepartitionFact"/>
      <Parameter name="CoarseSolver"                        type="string"   value="DirectSolver"/>
    </ParameterList>

  </ParameterList>
</ParameterList>

It is stored in ../../../test/tutorial/s5a.xml. In this example we define a smoothed aggregation transfer operator strategy (using the PgPFactory) for non-symmetric systems. The level smoother is chosen to be an over-relaxed symmetric Gauss–Seidel method. A direct solver is applied on the coarsest level. Please compare the building blocks in the xml file with Figure ref{fig:rebalanceddesignpgamg}. Be aware that the Nullspace variable now is also generated by the myRebalanceProlongatorFact.

Exercise 1 Choose option 1 in the problem menu of hands-on.py to run the Laplace 2D example on a \(300\times 300\) mesh. Change the solver to ../../../test/tutorial/s5a.xml. Use a reasonable number of processors. For demonstration purposes 4 processors should be fine for the \(300\times 300\) mesh. Run the example and check the screen output to see the effect of rebalancing. Try to visualize the ownership of the aggregates.

Note

The XML parameters in ../../../test/tutorial/s5a.xml write out the aggregation data for debugging. See the next tutorial for some more background information on aggregation and debugging.

You should observe a multigrid hierarchy as follows

Warning

Insert missing Screen output

printScreenOutput{s5a.txt_3.fragment_3.fragment}