#ifndef NGCORE_MPIWRAPPER_HPP #define NGCORE_MPIWRAPPER_HPP #include #include #include "array.hpp" #include "table.hpp" #include "exception.hpp" #include "profiler.hpp" #include "ngstream.hpp" #include "ng_mpi.hpp" namespace ngcore { #ifdef PARALLEL template struct MPI_typetrait { }; template <> struct MPI_typetrait { static NG_MPI_Datatype MPIType () { return NG_MPI_INT; } }; template <> struct MPI_typetrait { static NG_MPI_Datatype MPIType () { return NG_MPI_SHORT; } }; template <> struct MPI_typetrait { static NG_MPI_Datatype MPIType () { return NG_MPI_CHAR; } }; template <> struct MPI_typetrait { static NG_MPI_Datatype MPIType () { return NG_MPI_CHAR; } }; template <> struct MPI_typetrait { static NG_MPI_Datatype MPIType () { return NG_MPI_CHAR; } }; template <> struct MPI_typetrait { static NG_MPI_Datatype MPIType () { return NG_MPI_UINT64_T; } }; template <> struct MPI_typetrait { static NG_MPI_Datatype MPIType () { return NG_MPI_DOUBLE; } }; template <> struct MPI_typetrait> { static NG_MPI_Datatype MPIType () { return NG_MPI_CXX_DOUBLE_COMPLEX; } }; template <> struct MPI_typetrait { static NG_MPI_Datatype MPIType () { return NG_MPI_C_BOOL; } }; template struct MPI_typetrait> { static NG_MPI_Datatype MPIType () { static NG_MPI_Datatype NG_MPI_T = 0; if (!NG_MPI_T) { NG_MPI_Type_contiguous ( S, MPI_typetrait::MPIType(), &NG_MPI_T); NG_MPI_Type_commit ( &NG_MPI_T ); } return NG_MPI_T; } }; template ::MPIType())> inline NG_MPI_Datatype GetMPIType () { return MPI_typetrait::MPIType(); } template inline NG_MPI_Datatype GetMPIType (T &) { return GetMPIType(); } inline void MyMPI_WaitAll (FlatArray requests) { static Timer t("MPI - WaitAll"); RegionTimer reg(t); if (!requests.Size()) return; NG_MPI_Waitall (requests.Size(), requests.Data(), NG_MPI_STATUSES_IGNORE); } inline int MyMPI_WaitAny (FlatArray requests) { int nr; NG_MPI_Waitany (requests.Size(), requests.Data(), &nr, NG_MPI_STATUS_IGNORE); return nr; } class NgMPI_Comm { protected: NG_MPI_Comm comm; bool valid_comm; int * refcount; int rank, size; public: NgMPI_Comm () : valid_comm(false), refcount(nullptr), rank(0), size(1) { ; } NgMPI_Comm (NG_MPI_Comm _comm, bool owns = false) : comm(_comm), valid_comm(true) { int flag; NG_MPI_Initialized (&flag); if (!flag) { valid_comm = false; refcount = nullptr; rank = 0; size = 1; return; } if (!owns) refcount = nullptr; else refcount = new int{1}; NG_MPI_Comm_rank(comm, &rank); NG_MPI_Comm_size(comm, &size); } NgMPI_Comm (const NgMPI_Comm & c) : comm(c.comm), valid_comm(c.valid_comm), refcount(c.refcount), rank(c.rank), size(c.size) { if (refcount) (*refcount)++; } NgMPI_Comm (NgMPI_Comm && c) : comm(c.comm), valid_comm(c.valid_comm), refcount(c.refcount), rank(c.rank), size(c.size) { c.refcount = nullptr; } ~NgMPI_Comm() { if (refcount) if (--(*refcount) == 0) NG_MPI_Comm_free(&comm); } bool ValidCommunicator() const { return valid_comm; } NgMPI_Comm & operator= (const NgMPI_Comm & c) { if (refcount) if (--(*refcount) == 0) NG_MPI_Comm_free(&comm); refcount = c.refcount; if (refcount) (*refcount)++; comm = c.comm; valid_comm = c.valid_comm; size = c.size; rank = c.rank; return *this; } class InvalidCommException : public Exception { public: InvalidCommException() : Exception("Do not have a valid communicator") { ; } }; operator NG_MPI_Comm() const { if (!valid_comm) throw InvalidCommException(); return comm; } int Rank() const { return rank; } int Size() const { return size; } void Barrier() const { static Timer t("MPI - Barrier"); RegionTimer reg(t); if (size > 1) NG_MPI_Barrier (comm); } /** --- blocking P2P --- **/ template())> void Send (T & val, int dest, int tag) const { NG_MPI_Send (&val, 1, GetMPIType(), dest, tag, comm); } void Send (const std::string & s, int dest, int tag) const { NG_MPI_Send( const_cast (&s[0]), s.length(), NG_MPI_CHAR, dest, tag, comm); } template())> void Send(FlatArray s, int dest, int tag) const { NG_MPI_Send (s.Data(), s.Size(), GetMPIType(), dest, tag, comm); } template())> void Recv (T & val, int src, int tag) const { NG_MPI_Recv (&val, 1, GetMPIType(), src, tag, comm, NG_MPI_STATUS_IGNORE); } void Recv (std::string & s, int src, int tag) const { NG_MPI_Status status; int len; NG_MPI_Probe (src, tag, comm, &status); NG_MPI_Get_count (&status, NG_MPI_CHAR, &len); // s.assign (len, ' '); s.resize (len); NG_MPI_Recv( &s[0], len, NG_MPI_CHAR, src, tag, comm, NG_MPI_STATUS_IGNORE); } template ())> void Recv (FlatArray s, int src, int tag) const { NG_MPI_Recv (s.Data(), s.Size(), GetMPIType (), src, tag, comm, NG_MPI_STATUS_IGNORE); } template ())> void Recv (Array & s, int src, int tag) const { NG_MPI_Status status; int len; const NG_MPI_Datatype NG_MPI_T = GetMPIType (); NG_MPI_Probe (src, tag, comm, &status); NG_MPI_Get_count (&status, NG_MPI_T, &len); s.SetSize (len); NG_MPI_Recv (s.Data(), len, NG_MPI_T, src, tag, comm, NG_MPI_STATUS_IGNORE); } /** --- non-blocking P2P --- **/ template())> NG_MPI_Request ISend (T & val, int dest, int tag) const { NG_MPI_Request request; NG_MPI_Isend (&val, 1, GetMPIType(), dest, tag, comm, &request); return request; } template())> NG_MPI_Request ISend (FlatArray s, int dest, int tag) const { NG_MPI_Request request; NG_MPI_Isend (s.Data(), s.Size(), GetMPIType(), dest, tag, comm, &request); return request; } template())> NG_MPI_Request IRecv (T & val, int dest, int tag) const { NG_MPI_Request request; NG_MPI_Irecv (&val, 1, GetMPIType(), dest, tag, comm, &request); return request; } template())> NG_MPI_Request IRecv (FlatArray s, int src, int tag) const { NG_MPI_Request request; NG_MPI_Irecv (s.Data(), s.Size(), GetMPIType(), src, tag, comm, &request); return request; } /** --- collectives --- **/ template ())> T Reduce (T d, const NG_MPI_Op & op, int root = 0) const { static Timer t("MPI - Reduce"); RegionTimer reg(t); if (size == 1) return d; T global_d; NG_MPI_Reduce (&d, &global_d, 1, GetMPIType(), op, root, comm); return global_d; } template ())> T AllReduce (T d, const NG_MPI_Op & op) const { static Timer t("MPI - AllReduce"); RegionTimer reg(t); if (size == 1) return d; T global_d; NG_MPI_Allreduce ( &d, &global_d, 1, GetMPIType(), op, comm); return global_d; } template ())> void AllReduce (FlatArray d, const NG_MPI_Op & op) const { static Timer t("MPI - AllReduce Array"); RegionTimer reg(t); if (size == 1) return; NG_MPI_Allreduce (NG_MPI_IN_PLACE, d.Data(), d.Size(), GetMPIType(), op, comm); } template ())> void Bcast (T & s, int root = 0) const { if (size == 1) return; static Timer t("MPI - Bcast"); RegionTimer reg(t); NG_MPI_Bcast (&s, 1, GetMPIType(), root, comm); } template void Bcast (Array & d, int root = 0) { if (size == 1) return; int ds = d.Size(); Bcast (ds, root); if (Rank() != root) d.SetSize (ds); if (ds != 0) NG_MPI_Bcast (d.Data(), ds, GetMPIType(), root, comm); } void Bcast (std::string & s, int root = 0) const { if (size == 1) return; int len = s.length(); Bcast (len, root); if (rank != 0) s.resize (len); NG_MPI_Bcast (&s[0], len, NG_MPI_CHAR, root, comm); } template void AllToAll (FlatArray send, FlatArray recv) const { NG_MPI_Alltoall (send.Data(), 1, GetMPIType(), recv.Data(), 1, GetMPIType(), comm); } template void ScatterRoot (FlatArray send) const { if (size == 1) return; NG_MPI_Scatter (send.Data(), 1, GetMPIType(), NG_MPI_IN_PLACE, -1, GetMPIType(), 0, comm); } template void Scatter (T & recv) const { if (size == 1) return; NG_MPI_Scatter (NULL, 0, GetMPIType(), &recv, 1, GetMPIType(), 0, comm); } template void GatherRoot (FlatArray recv) const { recv[0] = T(0); if (size == 1) return; NG_MPI_Gather (NG_MPI_IN_PLACE, 1, GetMPIType(), recv.Data(), 1, GetMPIType(), 0, comm); } template void Gather (T send) const { if (size == 1) return; NG_MPI_Gather (&send, 1, GetMPIType(), NULL, 1, GetMPIType(), 0, comm); } template void AllGather (T val, FlatArray recv) const { if (size == 1) { recv[0] = val; return; } NG_MPI_Allgather (&val, 1, GetMPIType(), recv.Data(), 1, GetMPIType(), comm); } template void ExchangeTable (DynamicTable & send_data, DynamicTable & recv_data, int tag) { Array send_sizes(size); Array recv_sizes(size); for (int i = 0; i < size; i++) send_sizes[i] = send_data[i].Size(); AllToAll (send_sizes, recv_sizes); recv_data = DynamicTable (recv_sizes, true); Array requests; for (int dest = 0; dest < size; dest++) if (dest != rank && send_data[dest].Size()) requests.Append (ISend (FlatArray(send_data[dest]), dest, tag)); for (int dest = 0; dest < size; dest++) if (dest != rank && recv_data[dest].Size()) requests.Append (IRecv (FlatArray(recv_data[dest]), dest, tag)); MyMPI_WaitAll (requests); } NgMPI_Comm SubCommunicator (FlatArray procs) const { NG_MPI_Comm subcomm; NG_MPI_Group gcomm, gsubcomm; NG_MPI_Comm_group(comm, &gcomm); NG_MPI_Group_incl(gcomm, procs.Size(), procs.Data(), &gsubcomm); NG_MPI_Comm_create_group(comm, gsubcomm, 4242, &subcomm); return NgMPI_Comm(subcomm, true); } }; // class NgMPI_Comm class MyMPI { bool initialized_by_me; public: MyMPI(int argc, char ** argv) { int is_init = -1; NG_MPI_Initialized(&is_init); if (!is_init) { NG_MPI_Init (&argc, &argv); initialized_by_me = true; } else initialized_by_me = false; NgMPI_Comm comm(NG_MPI_COMM_WORLD); NGSOStream::SetGlobalActive (comm.Rank() == 0); if (comm.Size() > 1) TaskManager::SetNumThreads (1); } ~MyMPI() { if (initialized_by_me) NG_MPI_Finalize (); } }; #else // PARALLEL class NG_MPI_Comm { int nr; public: NG_MPI_Comm (int _nr = 0) : nr(_nr) { ; } operator int() const { return nr; } bool operator== (NG_MPI_Comm c2) const { return nr == c2.nr; } }; static NG_MPI_Comm NG_MPI_COMM_WORLD = 12345, NG_MPI_COMM_NULL = 10000; typedef int NG_MPI_Op; typedef int NG_MPI_Datatype; typedef int NG_MPI_Request; enum { NG_MPI_SUM = 0, NG_MPI_MIN = 1, NG_MPI_MAX = 2, NG_MPI_LOR = 4711 }; inline void NG_MPI_Type_contiguous ( int, NG_MPI_Datatype, NG_MPI_Datatype*) { ; } inline void NG_MPI_Type_commit ( NG_MPI_Datatype * ) { ; } template struct MPI_typetrait { static NG_MPI_Datatype MPIType () { return -1; } }; template inline NG_MPI_Datatype GetMPIType () { return -1; } class NgMPI_Comm { public: NgMPI_Comm () { ; } NgMPI_Comm (NG_MPI_Comm _comm, bool owns = false) { ; } size_t Rank() const { return 0; } size_t Size() const { return 1; } bool ValidCommunicator() const { return false; } void Barrier() const { ; } operator NG_MPI_Comm() const { return NG_MPI_Comm(); } template void Send( T & val, int dest, int tag) const { ; } template void Send(FlatArray s, int dest, int tag) const { ; } template void Recv (T & val, int src, int tag) const { ; } template void Recv (FlatArray s, int src, int tag) const { ; } template void Recv (Array & s, int src, int tag) const { ; } template NG_MPI_Request ISend (T & val, int dest, int tag) const { return 0; } template NG_MPI_Request ISend (FlatArray s, int dest, int tag) const { return 0; } template NG_MPI_Request IRecv (T & val, int dest, int tag) const { return 0; } template NG_MPI_Request IRecv (FlatArray s, int src, int tag) const { return 0; } template T Reduce (T d, const NG_MPI_Op & op, int root = 0) const { return d; } template T AllReduce (T d, const NG_MPI_Op & op) const { return d; } template void AllReduce (FlatArray d, const NG_MPI_Op & op) const { ; } template void Bcast (T & s, int root = 0) const { ; } template void Bcast (Array & d, int root = 0) { ; } template void AllGather (T val, FlatArray recv) const { recv[0] = val; } template void ExchangeTable (DynamicTable & send_data, DynamicTable & recv_data, int tag) { ; } NgMPI_Comm SubCommunicator (FlatArray procs) const { return *this; } }; inline void MyMPI_WaitAll (FlatArray requests) { ; } inline int MyMPI_WaitAny (FlatArray requests) { return 0; } class MyMPI { public: MyMPI(int argc, char ** argv) { ; } }; #endif // PARALLEL } // namespace ngcore #endif // NGCORE_MPIWRAPPER_HPP