Communicator.hpp 6.78 KB
Newer Older
Praetorius, Simon's avatar
Praetorius, Simon committed
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250
#pragma once

#if HAVE_MPI
#include <array>
#include <iostream>
#include <list>
#include <map>
#include <memory>
#include <string>
#include <type_traits>
#include <vector>

#include <mpi.h>

#include <dune/functions/common/functionconcepts.hh>

#include <amdis/common/parallel/Request.hpp>
#include <amdis/common/parallel/RecvDynamicSize.hpp>
#include <amdis/common/parallel/MpiTraits.hpp>

namespace AMDiS { namespace Mpi
{
  class Communicator
  {
  public:

    /// Constructor, stores an MPI communicator, e.g. MPI_COMM_WORLD
    Communicator(MPI_Comm comm)
      : comm_(comm)
      , buffer_(1024)
    {
      MPI_Comm_size(comm_, &size_);
      MPI_Comm_rank(comm_, &rank_);

      MPI_Buffer_attach(buffer_.data(), sizeof(char)*buffer_.size());
    }


  public:

    operator MPI_Comm() const { return comm_; }

    int size() const { return size_; }
    int rank() const { return rank_; }

  public:

    // send mpi datatype
    template <class Data>
    void send(Data const& data, int to, int tag = 0) const;

    // send array of mpi datatypes
    template <class T>
    void send(T const* data, std::size_t size, int to, int tag = 0) const;

    template <class T, std::size_t N>
    void send(T const (&data)[N], int to, int tag = 0) const
    {
      send(&data[0], N, to, tag);
    }

    template <class T, std::size_t N>
    void send(std::array<T,N> const& array, int to, int tag = 0) const
    {
      send(array.data(), N, to, tag);
    }

    template <class T>
    void send(std::vector<T> const& vec, int to, int tag = 0) const;

    void send(std::string const& str, int to, int tag = 0) const
    {
      MPI_Send(to_void_ptr(str.data()), int(str.size()), MPI_CHAR, to, tag, comm_);
    }


    // -------------------------------------------------------------------------------------


    // send mpi datatype (non-blocking)
    template <class Data>
    Request isend(Data const& data, int to, int tag = 0) const;

    // send mpi datatype (non-blocking, buffered)
    template <class Data>
    Request ibsend(Data const& data, int to, int tag = 0) const;

    // send array of mpi datatypes (non-blocking)
    template <class Data>
    Request isend(Data const* data, std::size_t size, int to, int tag = 0) const;

    // send array of mpi datatypes (non-blocking, buffered)
    template <class Data>
    Request ibsend(Data const* data, std::size_t size, int to, int tag = 0) const;

    template <class T>
    Request isend(std::vector<T> const& vec, int to, int tag = 0) const;

    Request isend(std::string const& str, int to, int tag = 0) const
    {
      MPI_Request request;
      MPI_Isend(to_void_ptr(str.data()), int(str.size()), MPI_CHAR, to, tag, comm_, &request);
      return {request};
    }

    // -------------------------------------------------------------------------------------

    // receive mpi datatype
    template <class Data>
    MPI_Status recv(Data& data, int from, int tag = 0) const;

    // receive array of mpi datatypes
    template <class T>
    MPI_Status recv(T* data, std::size_t size, int from, int tag = 0) const;

    template <class T, std::size_t N>
    MPI_Status recv(T (&data)[N], int from, int tag = 0) const
    {
      return recv(data, N, from, tag);
    }

    template <class T, std::size_t N>
    MPI_Status recv(std::array<T,N>& data, int from, int tag = 0) const
    {
      return recv(data.data(), N, from, tag);
    }

    template <class T>
    MPI_Status recv(std::vector<T>& data, int from, int tag = 0) const;

    MPI_Status recv(std::string& str, int from, int tag = 0) const
    {
      MPI_Status status;
      MPI_Probe(from, tag, comm_, &status);

      int size = 0;
      MPI_Get_count(&status, MPI_CHAR, &size);

      str.resize(size);
      MPI_Recv(&str[0], size, MPI_CHAR, from, tag, comm_, MPI_STATUS_IGNORE);
      return status;
    }

    // -------------------------------------------------------------------------------------

    // receive mpi datatype
    template <class Data>
    Request irecv(Data& data, int from, int tag = 0) const;

    // receive array of mpi datatypes
    template <class Data>
    Request irecv(Data* data, std::size_t size, int from, int tag = 0) const;

    template <class T, std::size_t N>
    Request irecv(T (&data)[N], int from, int tag = 0) const
    {
      return irecv(&data[0], N, from, tag);
    }

    template <class T, std::size_t N>
    Request irecv(std::array<T,N>& data, int from, int tag = 0) const
    {
      return irecv(data.data(), N, from, tag);
    }

    template <class Receiver>
      std::enable_if_t<Dune::Functions::Concept::isCallable<Receiver,MPI_Status>(), Request>
    irecv(Receiver&& recv, int from, int tag = 0) const
    {
      return Request{ RecvDynamicSize(from, tag, comm_, std::forward<Receiver>(recv)) };
    }

    // receive vector of mpi datatypes
    // 1. until message received, call MPI_Iprobe to retrieve status and size of message
    // 2. resize data-vector
    // 3. receive data into vector
    template <class T>
    Request irecv(std::vector<T>& vec, int from, int tag = 0) const;

    Request irecv(std::string& str, int from, int tag = 0) const
    {
      return Request{RecvDynamicSize(from,tag,comm_,
        [comm=comm_,&str](MPI_Status status) -> MPI_Request
        {
          int size = 0;
          MPI_Get_count(&status, MPI_CHAR, &size);

          str.resize(size);
          MPI_Request req;
          MPI_Irecv(&str[0], size, MPI_CHAR, status.MPI_SOURCE, status.MPI_TAG, comm, &req);

          return req;
        }) };
    }


  protected:

    // free unused buffers
    void check_buffers() const
    {
      using Buffers = std::decay_t<decltype(buffers_)>;

      std::list<typename Buffers::iterator> remove;
      for (auto it = buffers_.begin(); it != buffers_.end(); ++it) {
        int flag;
        MPI_Request_get_status(it->first, &flag, MPI_STATUS_IGNORE);
        if (flag != 0)
          remove.push_back(it);
      }

      for (auto it : remove)
        buffers_.erase(it);
    }


    std::pair<MPI_Request, std::string>& make_buffer(MPI_Status status, std::size_t len) const
    {
      auto it = buffers_.emplace(buffers_.end(), MPI_Request{}, std::string(len,' '));
      buffers_iterators_[{status.MPI_SOURCE, status.MPI_TAG}] = it;
      return buffers_.back();
    }

    std::pair<MPI_Request, std::string>& get_buffer(MPI_Status status) const
    {
      auto it = buffers_iterators_[{status.MPI_SOURCE, status.MPI_TAG}];
      return *it;
    }

  protected:

    MPI_Comm comm_;

    int rank_ = 0;
    int size_ = 1;

    std::vector<char> buffer_;

    using BufferList = std::list< std::pair<MPI_Request, std::string> >;
    mutable BufferList buffers_;

    using BufferIter = BufferList::iterator;
    mutable std::map<std::pair<int,int>, BufferIter> buffers_iterators_;
  };

}} // end namespace AMDiS::Mpi

#endif // HAVE_MPI

#include <amdis/common/parallel/Communicator.inc.hpp>