Master/Worker-Pattern mit MPI

Das Master/Worker-Pattern wird bei MPI gerne eingesetzt, da bei vielen Problemen sich nicht ohne weiteres eine Arbeitsaufteilung durchführen lässt, die alle beteiligten Prozesse gleichmäßig auslastet. Stattdessen ist es sinnvoll, Arbeitsaufträge häppchenweise auszuliefern.

Hierbei verwaltet der Master den gesamten Arbeitsumfang, schickt zu Beginn die ersten kleinen Aufträge heraus und wartet dann auf Antworten. Sobald eine Antwort eintrifft, kann in der Antwort der nächste Auftrag vergeben werden. Wenn alles abgearbeitet ist, muss der Master dies den Workern mitteilen, worauf diese dann ihre Aktivität einstellen. Die Unterscheidung zwischen Aufträgen und Ende-Mitteilungen kann durch unterschiedliche Tags erfolgen. Bei einem Ende-Paket kann zusätzlich auch die Größe des zu übertragenden Arrays auf 0 gesetzt werden.

Aufgabe

Zu den offenen Fragen der Zahlentheorie gehört, ob bestimmte aufeinanderfolgende Konstellationen von Primzahlen unendlich oft auftreten und wenn ja, wie sie verteilt sind. Diese Fragen werden gerne empirisch untersucht und lassen sich hervorragend parallelisieren.

In der Vorlage gegeben ist eine MPI-Anwendung, die ein Intervall der natürlichen Zahlen \([N_1, N_2]\) und eine vorgegebene Primzahlkonstellation \(\left\{n_i\right\}_{i=1}^{k-1}\) mit \(n_i < n_{i+1}\) der Länge \(k > 1\) erhält. Gesucht und zurückzuliefern sind dann alle Primzahlkonstellationen \((p, p+n_1, \dots p+n_{k-1})\) mit \(N_1 \le p \le N_2\). Primzahlpaare werden beispielsweise mit \(\left\{2\right\}\) gesucht, Primzahlvierlinge mit \(\left\{2, 6, 8\right\}\).

Die Vorlage verteilt das Gesamtintervall gleichmäßig auf alle Worker. Entwickeln Sie diese so weiter, dass gemäß dem Master/Worker-Pattern kleinere Teilintervalle an die einzelnen Worker vergeben werden, so dass bis zum Ende eine gleichmäßige Auslastung sichergestellt ist.

Da gerade auch größere Primzahlen reizvoll sind, arbeitet die Vorlage bereits mit beliebig großen ganzen Zahlen aus der GMP-Bibliothek. In der Vorlesungsbibliothek unter /home/numerik/pub/hpc/session20 gibt es in hpc/gmp/integer.h mit der Klasse ExportedInteger etwas Unterstützung, um diese Zahlen in int-Arrays zu verwandeln, die dann ganz einfach in MPI übertragen werden können. In der Vorlage erledigen das die Funktionen send_integer und receive_integer. Die probalistische Primzahlsuche selbst findet sich in primes.hpp und primes.cpp.

Beim Zusammenbau des Programms sind noch die Bibliotheksoptionen -lgmpxx -lgmp am Ende der Kommandozeile mit anzugeben.

#ifndef PRIMES_H
#define PRIMES_H

#include <hpc/gmp/integer.h>

bool search_prime_constellation(
         hpc::gmp::Integer& start_interval,
         hpc::gmp::Integer& end_interval,
         unsigned int k,           // size of prime constellation
         unsigned int offsets[],   // k-1 offsets
         hpc::gmp::Integer& first_prime);    // found result (if returning true)

#endif
#include "primes.hpp"

using namespace hpc::gmp;

bool search_prime_constellation(
         Integer& start_interval,
         Integer& end_interval,
         unsigned int k,         // size of prime constellation
         unsigned int offsets[], // k-1 offsets
         Integer& first_prime) {   // found result (if returning true)
   Integer p = start_interval;
   Integer q;
   bool found = false;
   for(;;) {
      // lookout for the next prime in the interval
      mpz_nextprime(p.get_mpz_t(), p.get_mpz_t());
      if (p > end_interval) break;
      // p is apparently prime, check for an constellation
      found = true;
      for (int i = 0; i < k-1; ++i) {
         unsigned int offset = offsets[i];
         q = p + offset;
         if (mpz_probab_prime_p(q.get_mpz_t(), 10) < 1) {
            found = false;
            break; // not a prime
         }
      }
      if (found) break;
   }
   if (found) {
      first_prime = p;
   }
   return found;
}
#include <mpi.h>
#include <iostream>
#include <cassert>
#include <cstdlib>
#include "primes.hpp"
#include <hpc/gmp/integer.h>
#include <hpc/aux/slices.h>

using namespace std;
using namespace hpc::aux;
using namespace hpc::gmp;

void send_integer(const Integer& value, int dest, int tag) {
   assert(dest >= 0);
   ExportedInteger exp_value(value);
   int len = (int) exp_value.length();
   int header[2] = {exp_value.sgn(), len};
   MPI_Send(header, 2, MPI_INT, dest, tag, MPI_COMM_WORLD);
   MPI_Send(exp_value.words, len, MPI_INT, dest, tag, MPI_COMM_WORLD);
}

void send_finish(int dest) {
   MPI_Send(nullptr, 0, MPI_INT, dest, 1, MPI_COMM_WORLD);
}

bool receive_integer(Integer& value, int& source, int& tag) {
   int lenval;
   MPI_Status status;
   int header[2];
   MPI_Recv(header, 2, MPI_INT, source, tag, MPI_COMM_WORLD, &status);
   tag = status.MPI_TAG;
   source = status.MPI_SOURCE;
   if (tag) return false;
   int sgn = header[0]; unsigned int len = header[1];
   ExportedInteger exp_value(sgn, len);
   MPI_Recv(exp_value.words, len, MPI_INT, source, tag, MPI_COMM_WORLD, &status);
   value = exp_value.get();
   return true;
}

char* progname;

void usage() {
   cerr << "Usage: " << progname << " N1 N2 {n_i} " << endl;
   exit(1);
}

void primes_master(int nofworkers, int argc, char** argv) {
   progname = *argv++; --argc;
   if (argc < 3) usage();
   Integer start(*argv++); --argc;
   Integer end(*argv++); --argc;

   int k = argc + 1;
   unsigned int* offsets = new unsigned int[k-1];
   for (int i = 0; i < k-1; ++i) {
      offsets[i] = atoi(*argv++); --argc;
   }

   MPI_Bcast(&k, 1, MPI_UNSIGNED, 0, MPI_COMM_WORLD);
   MPI_Bcast(offsets, k-1, MPI_UNSIGNED, 0, MPI_COMM_WORLD);

   for (int worker = 1; worker <= nofworkers; ++worker) {
      send_integer(start, worker, 0);
      send_integer(end, worker, 0);
   }
   int running_workers = nofworkers;
   while (running_workers > 0) {
      int source = MPI_ANY_SOURCE; int tag = MPI_ANY_TAG;
      Integer prime;
      bool ok = receive_integer(prime, source, tag);
      if (!ok) {
         std::cout << source << " has finished" << std::endl;
         --running_workers; continue;
      }
      std::cout << source << ": " << prime << std::endl;
   }
}

void primes_worker(int nofworkers, int rank) {
   int k;
   MPI_Bcast(&k, 1, MPI_UNSIGNED, 0, MPI_COMM_WORLD);
   unsigned int* offsets = new unsigned int[k-1];
   MPI_Bcast(offsets, k-1, MPI_UNSIGNED, 0, MPI_COMM_WORLD);

   Integer global_start, global_end;
   int source = 0, tag = 0;
   receive_integer(global_start, source, tag);
   receive_integer(global_end, source, tag);

   Slices<Integer> slices(nofworkers, global_end - global_start);
   Integer start = slices.offset(rank);
   Integer end = start + slices.size(rank);

   Integer prime;
   while (search_prime_constellation(start, end, k, offsets, prime)) {
      send_integer(prime, 0, 0);
      start = prime;
      start += 1;
   }
   send_finish(0);
}

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

   int rank; MPI_Comm_rank(MPI_COMM_WORLD, &rank);
   int nofworkers; MPI_Comm_size(MPI_COMM_WORLD, &nofworkers);
   --nofworkers; assert(nofworkers > 0);

   if (rank == 0) {
      primes_master(nofworkers, argc, argv);
   } else {
      primes_worker(nofworkers, rank-1);
   }

   MPI_Finalize();
}

Am Ende können Sie auch gerne Ihre Lösung wieder zu einem tar-Archiv zusammenpacken und mit submit einreichen:

submit hpc session20 session20.tar