Scatter/Gather

Zu den häufigen Kommunikationsszenarien analog zum Fork-und-Join-Pattern gehören

Dies ließe sich mit MPI_Send und MPI_Recv erledigen, wobei ein ausgewählter Prozess (typischerweise mit rank 0) dann wiederholt MPI_Send bzw. MPI_Recv aufrufen müsste.

Die sich wegen der begrenzten Pufferung aufaddierenden Latenzzeiten lassen sich reduzieren, wenn das Parallelisierungspotential durch entsprechende Operationen ausgenutzt wird. Die Vorgehensweise ist ähnlich wie bei MPI_Bcast, nur dass hier für jeden Prozess andere Inhalte versandt werden.

Entsprechende Operationen stehen auch in MPI mit MPI_Scatter und MPI_Gather zur Verfügung:

int MPI_Scatter(const void* sendbuf, int sendcount, MPI_Datatype sendtype,
                void* recvbuf, int recvcount, MPI_Datatype recvtype,
                int root, MPI_Comm comm);

int MPI_Gather(const void* sendbuf, int sendcount, MPI_Datatype sendtype,
               void* recvbuf, int recvcount, MPI_Datatype recvtype,
               int root, MPI_Comm comm);

Bei beiden Operationen ist root der rank des Prozesses, der eine Matrix oder einen Vektor verteilt bzw. wieder einsammelt. Normalerweise ist das die 0.

Bei MPI_Scatter erhält jeder der Prozesse (einschließlich root!) genau sendcount Elemente. Genauer, der Prozess mit dem Rank i erhält die Elemente \(i \cdot sendcount, \dots, (i+1) \cdot sendcount - 1\). Normalerweise entspricht sendcount dem Wert von recvcount, aber prinzipiell darf recvcount größer sein. Wie üblich müssen sendtype und recvtype kompatibel sein.

Analog legt bei MPI_Gather der recvcount fest, wieviel Elemente der Prozess mit dem Rank root von jedem der Prozesse (einschließlich sich selbst) erhält.

Ähnlich wie bei MPI_Bcast kann der gleiche Aufruf sowohl von root als auch von den anderen ausgeführt werden. Wenn es jedoch ungeschickt erscheint, die große Matrix in allen Prozessen anzulegen, kann es auch getrennt werden. Bei MPI_Scatter werden von den Nicht-Root-Prozessen die Parameter sendbuf, sendcount und sendtype ignoriert. Analog werden bei MPI_Gather die Parameter recvbuf, revcount und recvtype nur ausgewertet, wenn der eigene Rank root ist. Wenn der Code für beide Fälle getrennt ist, kann entsprechend nullptr bzw. 0 für die ungenutzten Parameter verwendet werden.

Aufgabe

Folgendes Testprogramm soll eine in Prozess 0 initialisierte Matrix zeilenweise auf die anderen Prozesse verteilen. Dabei erhält jeder Prozess share Zeilen. Nachdem jeder Prozess die Matrix bearbeitet hat, sollen die veränderten Resultate eingesammelt und ausgegeben werden.

Vorlage

#include <mpi.h>
#include <hpc/mpi/matrix.h>
#include <hpc/matvec/gematrix.h>
#include <hpc/matvec/print.h>
#include <hpc/matvec/apply.h>

int main(int argc, char** argv) {
   MPI_Init(&argc, &argv);

   int nof_processes; MPI_Comm_size(MPI_COMM_WORLD, &nof_processes);
   int rank; MPI_Comm_rank(MPI_COMM_WORLD, &rank);

   using namespace hpc::matvec;
   using namespace hpc::mpi;

   using Matrix = GeMatrix<double>;
   using Index = typename Matrix::Index;
   int share = 3;
   int num_rows = nof_processes * share;
   int num_cols = 5;

   Matrix B(share, num_cols, StorageOrder::RowMajor); /* individual share */

   if (rank == 0) {
      Matrix A(num_rows, num_cols); /* entire matrix */
      apply(A, [](double& val, Index i, Index j) -> void {
         val = i * 100 + j;
      });

      /* using MPI_Scatter: scatter A / receive our share into B */

      apply(B, [=](double& val, Index i, Index j) -> void {
         val += 10000 * (rank + 1);
      });

      /* using MPI_Gather: gather into A / send our share from B */

      print(A, "A");
   } else {

      /* using MPI_Scatter: receive our share into B */

      apply(B, [=](double& val, Index i, Index j) -> void {
         val += 10000 * (rank + 1);
      });

      /* using MPI_Gather: send our share from B */
   }
   MPI_Finalize();
}