Übertragung eigener Datentypen mit MPI

Es lohnt sich, noch einmal einen Blick auf MPI_Send und MPI_Recv zu werfen:

int MPI_Send(void* buf, int count, MPI_Datatype datatype,
             int dest, int tag, MPI_Comm comm);
int MPI_Recv(void* buf, int count, MPI_Datatype datatype,
             int source, int tag, MPI_Comm comm,
             MPI_Status* status);

Die 1:1-Nachrichtenaustausch-Operationen MPI_Send und MPI_Recv erlauben zwar die Übertragung von Arrays, bestehen aber darauf, dass die einzelnen Array-Elemente unmittelbar hintereinander im Speicher liegen. Somit können wir diese Operationen nicht unmittelbar nutzen, um einen beliebigen DenseVectorView oder einen GeMatrixView zu übertragen, da die Abstände, konfiguierbar durch inc bzw. incRow und incCol, beliebig sein können.

MPI bietet aber vielfältige Möglichkeiten Datentypen selbst zu definieren und es bietet sich an, passende MPI_Datatype-Objekte für unsere Vektor- und Matrixtypen zu erzeugen. Es ist zuvor jedoch sinnvoll den Aufbau eines Datentyps in MPI anzusehen:

Es gibt die Menge der Basistypen \(BT\) in MPI, der beispielsweise MPI_DOUBLE oder MPI_INT angehören. Ein Datentyp \(T\) mit der Kardinalität \(n\) ist in der MPI-Bibliothek eine Sequenz von Tupeln \(\left\{ (bt_1, o_1), (bt_2, o_2), \dots, (bt_n, o_n) \right\}\), mit \(bt_i \in BT\) und den zugehörigen Offsets \(o_i \in \mathbb{N}_0\) für \(i = 1, \dots, n\). Die Offsets geben die relative Position der jeweiligen Basiskomponenten in Bytes zur Anfangsadresse an.

Bezüglich der Kompatibilität bei MPI_Send und MPI_Recv werden zwei Datentypen \(T_1\) und \(T_2\) genau dann als kompatibel erachtet, falls die beiden Kardinalitäten \(n_1\) und \(n_2\) gleich sind und \(bt_{1_i} = bt_{2_i}\) für alle \(i = 1, \dots, n_1\) gilt. Bei \ident{MPI\Send} sind Überlappungen zulässig, bei \ident{MPI\Recv} haben sie einen undefinierten Effekt.

Einige Beispiele hierzu:

Da bezüglich der Kompatibilität die Offsets ignoriert werden, lässt sich bei der Übertragung problemlos ein Spalten- in einen Zeilenvektor konvertieren oder eine Matrix transponieren.

Für die Konstruktion von Typen gibt es zahlreiche Funktionen in MPI, die letztlich immer die beschriebenen Tupel-Sequenzen erzeugen. Ein erstes Beispiel hierzu sei MPI_Type_vector:

int MPI_Type_vector(int count, int blocklength, int stride,
                    MPI_Datatype oldtype, MPI_Datatype* newtype);

Der Elementtyp wird hier mit oldtype spezifiziert. Die blocklength spezifiziert, wieviel Elemente dieses Typs immer unmittelbar hintereinander liegen. Diese bilden einen Block. Der Offset zweier aufeinanderfolgenden Blöcke wird mit stride spezifiziert -- dieser Wert wird implizit mit der Größe des Elementtyps multipliziert. Insgesamt umfasst der Datentyp count Blöcke.

Bevor ein neu erzeugter Datentyp bei einer Übertragung wie beispielsweise mit MPI_Send und MPI_Recv verwendet werden kann, muss der Datentyp intern in eine flache Repräsentierung entsprechend der Tupel-Sequenz konvertiert werden. Dies geschieht mit MPI_Type_commit:

int MPI_Type_commit(MPI_Datatype* datatype);

Hier wird ein Zeiger auf den Datentyp übergeben, d.h. datatype kann anschließend einen neuen Wert haben.

Aufgabe

In der folgenden Vorlage ist ein kleiner Test zum Übertragen von Vektoren vorbereitet. Erzeugen Sie jeweils die notwendigen MPI-Datentypen, so dass jeweils einer der verwendeten Vektoren mit einer einzigen MPI_Send- bzw. MPI_Recv-Operation übertragen werden kann.

Hierbei lohnt es sich, eine kleine Funktion zu schreiben, die den MPI_Datatype für einen Vektor bestimmt. So könnte eine entsprechende Template-Funktion aussehen:

template<typename Vector>
typename std::enable_if<hpc::matvec::IsDenseVector<Vector>::value,
   MPI_Datatype>::type
get_type(const Vector& vector) {
   using ElementType = typename Vector::ElementType;
   MPI_Datatype datatype;
   /* ... */
   return datatype;
}

Wenn Sie einen allgemeinen Weg wünschen, um dem Element-Typen wiederum den konkreten MPI-Datentyp zuzuordnen, empfiehlt sich die Verwendung von #include <hpc/mpi/fundamental.h>. Dort gibt es bereits eine get_type-Template-Funktion, die alle elementaren Datentypen von MPI unterstützt. So liefert hpc::mpi::get_type(x) für double x bereits MPI_DOUBLE.

Die Vorlesungsbibliothek findet sich unter /home/numerik/pub/hpc/session21 auf unseren Rechnern.

Vorlage

#include <cassert>
#include <cstdio>
#include <mpi.h>
#include <hpc/matvec/gematrix.h>
#include <hpc/matvec/densevector.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);
   assert(nof_processes == 2);

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

   unsigned int nof_rows = 3;
   unsigned int nof_cols = 7;

   if (rank == 0) {
      GeMatrix<double> A(nof_rows, nof_cols, StorageOrder::RowMajor);
      using Index = GeMatrix<double>::Index;
      apply(A, [](double& val, Index i, Index j) -> void {
         val = i * 100 + j;
      });
      auto row = A.row(2);
      auto col = A.col(0);

      MPI_Send(/* row to rank 1 */);
      MPI_Send(/* col to rank 1 */);

      /* receive it back for verification */
      DenseVector<double> vec1(nof_cols), vec2(nof_rows);
      MPI_Status status;
      MPI_Recv(/* vec1 from rank 1 */);
      MPI_Recv(/* vec2 from rank 1 */);

      /* verify it */
      apply(vec1, [=](double& val, Index i) {
         if (val != row(i)) {
            printf("verification failed for row(%d): %lg vs %lg\n",
               (int)i, val, row(i));
         }
      });
      apply(vec2, [=](double& val, Index i) {
         if (val != col(i)) {
            printf("verification failed for col(%d): %lg vs %lg\n",
               (int)i, val, col(i));
         }
      });
   } else {
      DenseVector<double> vec1(nof_cols), vec2(nof_rows);
      MPI_Status status;
      MPI_Recv(/* vec1 from rank 0 */);
      MPI_Recv(/* vec2 from rank 0 */);

      /* send it back for verification */
      MPI_Send(/* vec1 to rank 0 */);
      MPI_Send(/* vec2 to rank 0 */);
   }
   MPI_Finalize();
}