netgen/libsrc/core/mpi_wrapper.hpp

634 lines
17 KiB
C++
Raw Normal View History

2019-02-12 00:01:07 +05:00
#ifndef NGCORE_MPIWRAPPER_HPP
#define NGCORE_MPIWRAPPER_HPP
2023-09-06 10:29:16 +05:00
#include <array>
2024-05-13 16:43:53 +05:00
#include <complex>
2019-02-12 00:01:07 +05:00
2019-08-26 16:02:13 +05:00
#include "array.hpp"
2022-04-27 01:00:25 +05:00
#include "table.hpp"
2019-07-10 13:56:55 +05:00
#include "exception.hpp"
2020-08-04 19:29:59 +05:00
#include "profiler.hpp"
2022-04-27 01:00:25 +05:00
#include "ngstream.hpp"
2024-05-13 16:43:53 +05:00
#include "ng_mpi.hpp"
2019-02-12 00:01:07 +05:00
namespace ngcore
{
#ifdef PARALLEL
template <class T> struct MPI_typetrait { };
template <> struct MPI_typetrait<int> {
2024-05-13 16:43:53 +05:00
static NG_MPI_Datatype MPIType () { return NG_MPI_INT; } };
2019-02-12 00:01:07 +05:00
template <> struct MPI_typetrait<short> {
2024-05-13 16:43:53 +05:00
static NG_MPI_Datatype MPIType () { return NG_MPI_SHORT; } };
2019-02-12 00:01:07 +05:00
template <> struct MPI_typetrait<char> {
2024-05-13 16:43:53 +05:00
static NG_MPI_Datatype MPIType () { return NG_MPI_CHAR; } };
2019-02-12 00:01:07 +05:00
2020-08-19 17:50:11 +05:00
template <> struct MPI_typetrait<signed char> {
2024-05-13 16:43:53 +05:00
static NG_MPI_Datatype MPIType () { return NG_MPI_CHAR; } };
2020-08-19 17:50:11 +05:00
2019-02-12 01:37:00 +05:00
template <> struct MPI_typetrait<unsigned char> {
2024-05-13 16:43:53 +05:00
static NG_MPI_Datatype MPIType () { return NG_MPI_CHAR; } };
2019-02-12 01:37:00 +05:00
2019-02-12 00:01:07 +05:00
template <> struct MPI_typetrait<size_t> {
2024-05-13 16:43:53 +05:00
static NG_MPI_Datatype MPIType () { return NG_MPI_UINT64_T; } };
2019-02-12 00:01:07 +05:00
template <> struct MPI_typetrait<double> {
2024-05-13 16:43:53 +05:00
static NG_MPI_Datatype MPIType () { return NG_MPI_DOUBLE; } };
template <> struct MPI_typetrait<std::complex<double>> {
static NG_MPI_Datatype MPIType () { return NG_MPI_CXX_DOUBLE_COMPLEX; } };
2019-02-12 00:01:07 +05:00
template <> struct MPI_typetrait<bool> {
2024-05-13 16:43:53 +05:00
static NG_MPI_Datatype MPIType () { return NG_MPI_C_BOOL; } };
2019-02-12 00:01:07 +05:00
2023-09-06 00:53:55 +05:00
2023-09-07 23:21:44 +05:00
template<typename T, size_t S>
2023-09-06 00:53:55 +05:00
struct MPI_typetrait<std::array<T,S>>
{
2024-05-13 16:43:53 +05:00
static NG_MPI_Datatype MPIType ()
2023-09-06 00:53:55 +05:00
{
2024-05-13 16:43:53 +05:00
static NG_MPI_Datatype NG_MPI_T = 0;
if (!NG_MPI_T)
2023-09-06 00:53:55 +05:00
{
2024-05-13 16:43:53 +05:00
NG_MPI_Type_contiguous ( S, MPI_typetrait<T>::MPIType(), &NG_MPI_T);
NG_MPI_Type_commit ( &NG_MPI_T );
2023-09-06 00:53:55 +05:00
}
2024-05-13 16:43:53 +05:00
return NG_MPI_T;
2023-09-06 00:53:55 +05:00
}
};
2019-02-12 00:01:07 +05:00
template <class T, class T2 = decltype(MPI_typetrait<T>::MPIType())>
2024-05-13 16:43:53 +05:00
inline NG_MPI_Datatype GetMPIType () {
2019-02-12 00:01:07 +05:00
return MPI_typetrait<T>::MPIType();
}
2020-08-19 17:50:11 +05:00
template <class T>
2024-05-13 16:43:53 +05:00
inline NG_MPI_Datatype GetMPIType (T &) {
2020-08-19 17:50:11 +05:00
return GetMPIType<T>();
}
2019-02-12 00:01:07 +05:00
2024-11-28 01:16:36 +05:00
class NgMPI_Request
{
NG_MPI_Request request;
public:
NgMPI_Request (NG_MPI_Request requ) : request{requ} { }
NgMPI_Request (const NgMPI_Request&) = delete;
NgMPI_Request (NgMPI_Request&&) = default;
~NgMPI_Request () { NG_MPI_Wait (&request, NG_MPI_STATUS_IGNORE); }
void Wait() { NG_MPI_Wait (&request, NG_MPI_STATUS_IGNORE); }
operator NG_MPI_Request() &&
{
auto tmp = request;
request = NG_MPI_REQUEST_NULL;
return tmp;
}
};
class NgMPI_Requests
{
Array<NG_MPI_Request> requests;
public:
NgMPI_Requests() = default;
~NgMPI_Requests() { WaitAll(); }
NgMPI_Requests & operator+= (NgMPI_Request && r)
{
requests += NG_MPI_Request(std::move(r));
return *this;
}
NgMPI_Requests & operator+= (NG_MPI_Request r)
{
requests += r;
return *this;
}
void WaitAll()
{
static Timer t("NgMPI - WaitAll"); RegionTimer reg(t);
if (!requests.Size()) return;
NG_MPI_Waitall (requests.Size(), requests.Data(), NG_MPI_STATUSES_IGNORE);
}
int WaitAny ()
{
int nr;
NG_MPI_Waitany (requests.Size(), requests.Data(), &nr, NG_MPI_STATUS_IGNORE);
return nr;
}
};
2022-04-27 01:00:25 +05:00
2024-05-13 16:43:53 +05:00
inline void MyMPI_WaitAll (FlatArray<NG_MPI_Request> requests)
2022-04-27 01:00:25 +05:00
{
static Timer t("MPI - WaitAll"); RegionTimer reg(t);
if (!requests.Size()) return;
2024-05-13 16:43:53 +05:00
NG_MPI_Waitall (requests.Size(), requests.Data(), NG_MPI_STATUSES_IGNORE);
2022-04-27 01:00:25 +05:00
}
2024-05-13 16:43:53 +05:00
inline int MyMPI_WaitAny (FlatArray<NG_MPI_Request> requests)
2022-04-27 01:00:25 +05:00
{
int nr;
2024-05-13 16:43:53 +05:00
NG_MPI_Waitany (requests.Size(), requests.Data(), &nr, NG_MPI_STATUS_IGNORE);
2022-04-27 01:00:25 +05:00
return nr;
}
2019-02-12 00:01:07 +05:00
class NgMPI_Comm
{
2019-02-13 04:11:35 +05:00
protected:
2024-05-13 16:43:53 +05:00
NG_MPI_Comm comm;
2019-02-13 16:28:24 +05:00
bool valid_comm;
2019-02-12 00:01:07 +05:00
int * refcount;
2019-02-12 01:37:00 +05:00
int rank, size;
2019-02-12 00:01:07 +05:00
public:
2019-02-12 12:03:20 +05:00
NgMPI_Comm ()
2019-02-13 16:28:24 +05:00
: valid_comm(false), refcount(nullptr), rank(0), size(1)
2019-02-12 12:03:20 +05:00
{ ; }
2024-05-13 16:43:53 +05:00
NgMPI_Comm (NG_MPI_Comm _comm, bool owns = false)
2019-02-13 16:28:24 +05:00
: comm(_comm), valid_comm(true)
2019-02-12 00:01:07 +05:00
{
2020-09-16 02:15:50 +05:00
int flag;
2024-05-13 16:43:53 +05:00
NG_MPI_Initialized (&flag);
2020-09-16 02:15:50 +05:00
if (!flag)
{
valid_comm = false;
refcount = nullptr;
rank = 0;
size = 1;
return;
}
2019-02-12 00:01:07 +05:00
if (!owns)
refcount = nullptr;
else
refcount = new int{1};
2019-02-12 01:37:00 +05:00
2024-05-13 16:43:53 +05:00
NG_MPI_Comm_rank(comm, &rank);
NG_MPI_Comm_size(comm, &size);
2019-02-12 00:01:07 +05:00
}
NgMPI_Comm (const NgMPI_Comm & c)
2019-02-13 16:28:24 +05:00
: comm(c.comm), valid_comm(c.valid_comm), refcount(c.refcount),
rank(c.rank), size(c.size)
2019-02-12 00:01:07 +05:00
{
if (refcount) (*refcount)++;
}
NgMPI_Comm (NgMPI_Comm && c)
2019-02-13 16:28:24 +05:00
: comm(c.comm), valid_comm(c.valid_comm), refcount(c.refcount),
rank(c.rank), size(c.size)
2019-02-12 00:01:07 +05:00
{
c.refcount = nullptr;
}
~NgMPI_Comm()
{
if (refcount)
if (--(*refcount) == 0)
2024-05-13 16:43:53 +05:00
NG_MPI_Comm_free(&comm);
2019-02-12 00:01:07 +05:00
}
2019-02-12 05:21:56 +05:00
bool ValidCommunicator() const
{
return valid_comm;
}
2019-02-12 05:21:56 +05:00
NgMPI_Comm & operator= (const NgMPI_Comm & c)
{
if (refcount)
if (--(*refcount) == 0)
2024-05-13 16:43:53 +05:00
NG_MPI_Comm_free(&comm);
2019-02-12 05:21:56 +05:00
refcount = c.refcount;
if (refcount) (*refcount)++;
2019-02-12 16:41:08 +05:00
comm = c.comm;
2019-02-13 16:28:24 +05:00
valid_comm = c.valid_comm;
2019-02-12 05:21:56 +05:00
size = c.size;
rank = c.rank;
return *this;
}
2019-02-13 16:28:24 +05:00
class InvalidCommException : public Exception {
public:
InvalidCommException() : Exception("Do not have a valid communicator") { ; }
};
2019-02-12 05:21:56 +05:00
2024-05-13 16:43:53 +05:00
operator NG_MPI_Comm() const {
2019-02-13 16:28:24 +05:00
if (!valid_comm) throw InvalidCommException();
return comm;
}
2019-02-12 00:01:07 +05:00
2019-02-13 02:11:55 +05:00
int Rank() const { return rank; }
int Size() const { return size; }
void Barrier() const {
2020-08-04 19:29:59 +05:00
static Timer t("MPI - Barrier"); RegionTimer reg(t);
2024-05-13 16:43:53 +05:00
if (size > 1) NG_MPI_Barrier (comm);
2019-02-13 02:11:55 +05:00
}
2019-02-12 01:37:00 +05:00
2019-02-13 04:11:35 +05:00
/** --- blocking P2P --- **/
2019-02-12 01:37:00 +05:00
template<typename T, typename T2 = decltype(GetMPIType<T>())>
2019-02-13 02:11:55 +05:00
void Send (T & val, int dest, int tag) const {
2024-05-13 16:43:53 +05:00
NG_MPI_Send (&val, 1, GetMPIType<T>(), dest, tag, comm);
2019-02-12 01:37:00 +05:00
}
2020-08-05 04:11:26 +05:00
void Send (const std::string & s, int dest, int tag) const {
2024-05-13 16:43:53 +05:00
NG_MPI_Send( const_cast<char*> (&s[0]), s.length(), NG_MPI_CHAR, dest, tag, comm);
2020-08-05 04:11:26 +05:00
}
2019-02-12 01:37:00 +05:00
2020-08-19 17:50:11 +05:00
template<typename T, typename TI, typename T2 = decltype(GetMPIType<T>())>
void Send(FlatArray<T,TI> s, int dest, int tag) const {
2024-05-13 16:43:53 +05:00
NG_MPI_Send (s.Data(), s.Size(), GetMPIType<T>(), dest, tag, comm);
2019-08-26 16:02:13 +05:00
}
2019-02-12 01:37:00 +05:00
template<typename T, typename T2 = decltype(GetMPIType<T>())>
2019-02-13 04:11:35 +05:00
void Recv (T & val, int src, int tag) const {
2024-05-13 16:43:53 +05:00
NG_MPI_Recv (&val, 1, GetMPIType<T>(), src, tag, comm, NG_MPI_STATUS_IGNORE);
2019-02-12 01:37:00 +05:00
}
2020-08-05 04:11:26 +05:00
void Recv (std::string & s, int src, int tag) const {
2024-05-13 16:43:53 +05:00
NG_MPI_Status status;
2020-08-05 04:11:26 +05:00
int len;
2024-05-13 16:43:53 +05:00
NG_MPI_Probe (src, tag, comm, &status);
NG_MPI_Get_count (&status, NG_MPI_CHAR, &len);
2020-08-05 04:11:26 +05:00
// s.assign (len, ' ');
s.resize (len);
2024-05-13 16:43:53 +05:00
NG_MPI_Recv( &s[0], len, NG_MPI_CHAR, src, tag, comm, NG_MPI_STATUS_IGNORE);
2020-08-05 04:11:26 +05:00
}
2020-08-19 17:50:11 +05:00
template <typename T, typename TI, typename T2 = decltype(GetMPIType<T>())>
void Recv (FlatArray <T,TI> s, int src, int tag) const {
2024-05-13 16:43:53 +05:00
NG_MPI_Recv (s.Data(), s.Size(), GetMPIType<T> (), src, tag, comm, NG_MPI_STATUS_IGNORE);
2019-08-26 16:02:13 +05:00
}
2020-08-19 17:50:11 +05:00
template <typename T, typename TI, typename T2 = decltype(GetMPIType<T>())>
void Recv (Array <T,TI> & s, int src, int tag) const
2019-08-26 16:02:13 +05:00
{
2024-05-13 16:43:53 +05:00
NG_MPI_Status status;
2019-08-26 16:02:13 +05:00
int len;
2024-05-13 16:43:53 +05:00
const NG_MPI_Datatype NG_MPI_T = GetMPIType<T> ();
NG_MPI_Probe (src, tag, comm, &status);
NG_MPI_Get_count (&status, NG_MPI_T, &len);
2019-08-26 16:02:13 +05:00
s.SetSize (len);
2024-05-13 16:43:53 +05:00
NG_MPI_Recv (s.Data(), len, NG_MPI_T, src, tag, comm, NG_MPI_STATUS_IGNORE);
2019-08-26 16:02:13 +05:00
}
2019-02-13 02:11:55 +05:00
2019-02-13 04:11:35 +05:00
/** --- non-blocking P2P --- **/
template<typename T, typename T2 = decltype(GetMPIType<T>())>
2024-05-13 16:43:53 +05:00
NG_MPI_Request ISend (T & val, int dest, int tag) const
2019-02-13 04:11:35 +05:00
{
2024-05-13 16:43:53 +05:00
NG_MPI_Request request;
NG_MPI_Isend (&val, 1, GetMPIType<T>(), dest, tag, comm, &request);
2019-02-13 04:11:35 +05:00
return request;
}
2019-08-28 17:09:51 +05:00
template<typename T, typename T2 = decltype(GetMPIType<T>())>
2024-05-13 16:43:53 +05:00
NG_MPI_Request ISend (FlatArray<T> s, int dest, int tag) const
2019-08-28 17:09:51 +05:00
{
2024-05-13 16:43:53 +05:00
NG_MPI_Request request;
NG_MPI_Isend (s.Data(), s.Size(), GetMPIType<T>(), dest, tag, comm, &request);
2019-08-28 17:09:51 +05:00
return request;
}
2019-02-13 04:11:35 +05:00
template<typename T, typename T2 = decltype(GetMPIType<T>())>
2024-05-13 16:43:53 +05:00
NG_MPI_Request IRecv (T & val, int dest, int tag) const
2019-02-13 04:11:35 +05:00
{
2024-05-13 16:43:53 +05:00
NG_MPI_Request request;
NG_MPI_Irecv (&val, 1, GetMPIType<T>(), dest, tag, comm, &request);
2019-02-13 04:11:35 +05:00
return request;
}
2019-08-26 16:02:13 +05:00
template<typename T, typename T2 = decltype(GetMPIType<T>())>
2024-05-13 16:43:53 +05:00
NG_MPI_Request IRecv (FlatArray<T> s, int src, int tag) const
2019-08-26 16:02:13 +05:00
{
2024-05-13 16:43:53 +05:00
NG_MPI_Request request;
NG_MPI_Irecv (s.Data(), s.Size(), GetMPIType<T>(), src, tag, comm, &request);
2019-08-26 16:02:13 +05:00
return request;
}
2019-02-13 02:11:55 +05:00
/** --- collectives --- **/
2019-02-13 04:11:35 +05:00
template <typename T, typename T2 = decltype(GetMPIType<T>())>
2024-05-13 16:43:53 +05:00
T Reduce (T d, const NG_MPI_Op & op, int root = 0) const
2019-02-13 04:11:35 +05:00
{
2020-08-04 19:29:59 +05:00
static Timer t("MPI - Reduce"); RegionTimer reg(t);
2019-02-13 04:11:35 +05:00
if (size == 1) return d;
T global_d;
2024-05-13 16:43:53 +05:00
NG_MPI_Reduce (&d, &global_d, 1, GetMPIType<T>(), op, root, comm);
2019-02-13 04:11:35 +05:00
return global_d;
}
2019-02-13 02:11:55 +05:00
template <typename T, typename T2 = decltype(GetMPIType<T>())>
2024-05-13 16:43:53 +05:00
T AllReduce (T d, const NG_MPI_Op & op) const
2019-02-13 02:11:55 +05:00
{
2020-08-04 19:29:59 +05:00
static Timer t("MPI - AllReduce"); RegionTimer reg(t);
2019-02-13 02:11:55 +05:00
if (size == 1) return d;
T global_d;
2024-05-13 16:43:53 +05:00
NG_MPI_Allreduce ( &d, &global_d, 1, GetMPIType<T>(), op, comm);
2019-02-13 02:11:55 +05:00
return global_d;
}
2022-03-05 15:42:00 +05:00
template <typename T, typename T2 = decltype(GetMPIType<T>())>
2024-05-13 16:43:53 +05:00
void AllReduce (FlatArray<T> d, const NG_MPI_Op & op) const
2022-03-05 15:42:00 +05:00
{
static Timer t("MPI - AllReduce Array"); RegionTimer reg(t);
2022-03-05 16:21:28 +05:00
if (size == 1) return;
2022-03-05 15:42:00 +05:00
2024-05-13 16:43:53 +05:00
NG_MPI_Allreduce (NG_MPI_IN_PLACE, d.Data(), d.Size(), GetMPIType<T>(), op, comm);
2022-03-05 15:42:00 +05:00
}
2019-02-13 02:11:55 +05:00
template <typename T, typename T2 = decltype(GetMPIType<T>())>
void Bcast (T & s, int root = 0) const {
2022-04-27 01:00:25 +05:00
if (size == 1) return;
2021-06-08 22:08:07 +05:00
static Timer t("MPI - Bcast"); RegionTimer reg(t);
2024-05-13 16:43:53 +05:00
NG_MPI_Bcast (&s, 1, GetMPIType<T>(), root, comm);
2019-02-13 02:11:55 +05:00
}
2022-04-27 01:00:25 +05:00
2024-11-26 17:29:14 +05:00
template <class T, size_t S>
void Bcast (std::array<T,S> & d, int root = 0) const
{
if (size == 1) return;
if (S != 0)
NG_MPI_Bcast (&d[0], S, GetMPIType<T>(), root, comm);
}
2022-04-27 01:00:25 +05:00
template <class T>
2024-11-26 17:29:14 +05:00
void Bcast (Array<T> & d, int root = 0) const
2022-04-27 01:00:25 +05:00
{
if (size == 1) return;
int ds = d.Size();
Bcast (ds, root);
if (Rank() != root) d.SetSize (ds);
if (ds != 0)
2024-05-13 16:43:53 +05:00
NG_MPI_Bcast (d.Data(), ds, GetMPIType<T>(), root, comm);
2022-04-27 01:00:25 +05:00
}
2019-02-13 02:11:55 +05:00
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);
2024-05-13 16:43:53 +05:00
NG_MPI_Bcast (&s[0], len, NG_MPI_CHAR, root, comm);
2019-02-13 02:11:55 +05:00
}
2024-11-26 17:29:14 +05:00
template <class T, size_t S>
2024-11-28 01:16:36 +05:00
NgMPI_Request IBcast (std::array<T,S> & d, int root = 0) const
2024-11-26 17:29:14 +05:00
{
NG_MPI_Request request;
NG_MPI_Ibcast (&d[0], S, GetMPIType<T>(), root, comm, &request);
return request;
}
template <class T>
2024-11-28 01:16:36 +05:00
NgMPI_Request IBcast (FlatArray<T> d, int root = 0) const
2024-11-26 17:29:14 +05:00
{
NG_MPI_Request request;
int ds = d.Size();
NG_MPI_Ibcast (d.Data(), ds, GetMPIType<T>(), root, comm, &request);
return request;
}
2020-07-31 12:57:19 +05:00
template <typename T>
void AllToAll (FlatArray<T> send, FlatArray<T> recv) const
{
2024-05-13 16:43:53 +05:00
NG_MPI_Alltoall (send.Data(), 1, GetMPIType<T>(),
2020-07-31 12:57:19 +05:00
recv.Data(), 1, GetMPIType<T>(), comm);
}
template <typename T>
void ScatterRoot (FlatArray<T> send) const
{
if (size == 1) return;
2024-05-13 16:43:53 +05:00
NG_MPI_Scatter (send.Data(), 1, GetMPIType<T>(),
NG_MPI_IN_PLACE, -1, GetMPIType<T>(), 0, comm);
}
template <typename T>
void Scatter (T & recv) const
{
if (size == 1) return;
2024-05-13 16:43:53 +05:00
NG_MPI_Scatter (NULL, 0, GetMPIType<T>(),
&recv, 1, GetMPIType<T>(), 0, comm);
}
2022-03-18 12:20:20 +05:00
template <typename T>
void GatherRoot (FlatArray<T> recv) const
{
recv[0] = T(0);
if (size == 1) return;
2024-05-13 16:43:53 +05:00
NG_MPI_Gather (NG_MPI_IN_PLACE, 1, GetMPIType<T>(),
2022-03-18 12:20:20 +05:00
recv.Data(), 1, GetMPIType<T>(), 0, comm);
}
template <typename T>
void Gather (T send) const
{
if (size == 1) return;
2024-05-13 16:43:53 +05:00
NG_MPI_Gather (&send, 1, GetMPIType<T>(),
2022-03-18 12:20:20 +05:00
NULL, 1, GetMPIType<T>(), 0, comm);
}
2020-08-17 23:28:00 +05:00
template <typename T>
void AllGather (T val, FlatArray<T> recv) const
{
if (size == 1)
{
recv[0] = val;
return;
}
2024-05-13 16:43:53 +05:00
NG_MPI_Allgather (&val, 1, GetMPIType<T>(),
2020-08-17 23:28:00 +05:00
recv.Data(), 1, GetMPIType<T>(),
comm);
}
2020-07-31 12:57:19 +05:00
2022-04-27 01:00:25 +05:00
template <typename T>
void ExchangeTable (DynamicTable<T> & send_data,
DynamicTable<T> & recv_data, int tag)
{
Array<int> send_sizes(size);
Array<int> 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<T> (recv_sizes, true);
2024-05-13 16:43:53 +05:00
Array<NG_MPI_Request> requests;
2022-04-27 01:00:25 +05:00
for (int dest = 0; dest < size; dest++)
if (dest != rank && send_data[dest].Size())
requests.Append (ISend (FlatArray<T>(send_data[dest]), dest, tag));
for (int dest = 0; dest < size; dest++)
if (dest != rank && recv_data[dest].Size())
requests.Append (IRecv (FlatArray<T>(recv_data[dest]), dest, tag));
MyMPI_WaitAll (requests);
}
2020-07-31 12:57:19 +05:00
2019-08-28 17:59:32 +05:00
NgMPI_Comm SubCommunicator (FlatArray<int> procs) const
{
2024-05-13 16:43:53 +05:00
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);
2019-08-28 17:59:32 +05:00
return NgMPI_Comm(subcomm, true);
}
2019-08-26 16:02:13 +05:00
}; // class NgMPI_Comm
2019-07-10 13:56:55 +05:00
#else // PARALLEL
2024-05-13 16:43:53 +05:00
class NG_MPI_Comm {
2019-02-12 03:13:12 +05:00
int nr;
public:
2024-05-13 16:43:53 +05:00
NG_MPI_Comm (int _nr = 0) : nr(_nr) { ; }
2019-02-12 03:13:12 +05:00
operator int() const { return nr; }
2024-05-13 16:43:53 +05:00
bool operator== (NG_MPI_Comm c2) const { return nr == c2.nr; }
2019-02-12 03:13:12 +05:00
};
2024-05-13 16:43:53 +05:00
static NG_MPI_Comm NG_MPI_COMM_WORLD = 12345, NG_MPI_COMM_NULL = 10000;
2019-02-13 02:25:22 +05:00
2024-05-13 16:43:53 +05:00
typedef int NG_MPI_Op;
typedef int NG_MPI_Datatype;
typedef int NG_MPI_Request;
2019-02-13 04:11:35 +05:00
2024-05-13 16:43:53 +05:00
enum { NG_MPI_SUM = 0, NG_MPI_MIN = 1, NG_MPI_MAX = 2, NG_MPI_LOR = 4711 };
2024-05-13 16:43:53 +05:00
inline void NG_MPI_Type_contiguous ( int, NG_MPI_Datatype, NG_MPI_Datatype*) { ; }
inline void NG_MPI_Type_commit ( NG_MPI_Datatype * ) { ; }
2023-09-06 10:49:27 +05:00
2023-09-06 11:32:32 +05:00
template <class T> struct MPI_typetrait {
2024-05-13 16:43:53 +05:00
static NG_MPI_Datatype MPIType () { return -1; }
2023-09-06 11:32:32 +05:00
};
template <class T, class T2=void>
2024-05-13 16:43:53 +05:00
inline NG_MPI_Datatype GetMPIType () { return -1; }
2024-11-28 01:16:36 +05:00
class NgMPI_Request { };
class NgMPI_Requests
{
public:
NgMPI_Requests operator+= (NgMPI_Request &&) { ; }
NgMPI_Requests operator+= (NG_MPI_Request r) { ; }
void WaitAll() { ; }
int WaitAny() { return 0; }
};
2019-02-12 00:01:07 +05:00
class NgMPI_Comm
{
public:
2019-02-12 12:03:20 +05:00
NgMPI_Comm () { ; }
2024-05-13 16:43:53 +05:00
NgMPI_Comm (NG_MPI_Comm _comm, bool owns = false) { ; }
2019-02-12 00:01:07 +05:00
size_t Rank() const { return 0; }
size_t Size() const { return 1; }
2020-08-02 14:33:11 +05:00
bool ValidCommunicator() const { return false; }
2019-02-13 02:11:55 +05:00
void Barrier() const { ; }
2024-05-13 16:43:53 +05:00
operator NG_MPI_Comm() const { return NG_MPI_Comm(); }
2019-02-12 01:37:00 +05:00
template<typename T>
2019-02-13 02:11:55 +05:00
void Send( T & val, int dest, int tag) const { ; }
2019-02-12 01:37:00 +05:00
template<typename T>
2019-08-26 16:02:13 +05:00
void Send(FlatArray<T> s, int dest, int tag) const { ; }
template<typename T>
void Recv (T & val, int src, int tag) const { ; }
template <typename T>
void Recv (FlatArray <T> s, int src, int tag) const { ; }
template <typename T>
void Recv (Array <T> & s, int src, int tag) const { ; }
2019-02-13 02:11:55 +05:00
2019-02-13 04:11:35 +05:00
template<typename T>
2024-05-13 16:43:53 +05:00
NG_MPI_Request ISend (T & val, int dest, int tag) const { return 0; }
2019-02-13 04:11:35 +05:00
2019-08-26 16:02:13 +05:00
template<typename T>
2024-05-13 16:43:53 +05:00
NG_MPI_Request ISend (FlatArray<T> s, int dest, int tag) const { return 0; }
2019-08-26 16:02:13 +05:00
2019-02-13 04:11:35 +05:00
template<typename T>
2024-05-13 16:43:53 +05:00
NG_MPI_Request IRecv (T & val, int dest, int tag) const { return 0; }
2019-02-13 04:11:35 +05:00
2019-08-26 16:02:13 +05:00
template<typename T>
2024-05-13 16:43:53 +05:00
NG_MPI_Request IRecv (FlatArray<T> s, int src, int tag) const { return 0; }
2019-08-26 16:02:13 +05:00
2019-02-13 04:11:35 +05:00
template <typename T>
2024-05-13 16:43:53 +05:00
T Reduce (T d, const NG_MPI_Op & op, int root = 0) const { return d; }
2019-02-13 04:11:35 +05:00
2019-02-13 02:16:53 +05:00
template <typename T>
2024-05-13 16:43:53 +05:00
T AllReduce (T d, const NG_MPI_Op & op) const { return d; }
2019-02-13 02:11:55 +05:00
2022-03-05 15:42:00 +05:00
template <typename T>
2024-05-13 16:43:53 +05:00
void AllReduce (FlatArray<T> d, const NG_MPI_Op & op) const { ; }
2022-03-05 15:42:00 +05:00
2019-02-13 02:16:53 +05:00
template <typename T>
2019-02-13 02:21:44 +05:00
void Bcast (T & s, int root = 0) const { ; }
2019-08-26 16:02:13 +05:00
2024-11-26 17:29:14 +05:00
template <class T, size_t S>
void Bcast (std::array<T,S> & d, int root = 0) const {}
2022-04-27 01:00:25 +05:00
template <class T>
2024-11-26 17:29:14 +05:00
void Bcast (Array<T> & d, int root = 0) const { ; }
template <class T, size_t S>
NG_MPI_Request IBcast (std::array<T,S> & d, int root = 0) const { return 0; }
2022-04-27 01:00:25 +05:00
2024-11-26 17:29:14 +05:00
template <class T>
NG_MPI_Request IBcast (FlatArray<T> d, int root = 0) const { return 0; }
2022-05-06 19:39:06 +05:00
template <typename T>
void AllGather (T val, FlatArray<T> recv) const
{
recv[0] = val;
}
2022-04-27 01:00:25 +05:00
template <typename T>
void ExchangeTable (DynamicTable<T> & send_data,
DynamicTable<T> & recv_data, int tag) { ; }
2019-08-26 16:02:13 +05:00
NgMPI_Comm SubCommunicator (FlatArray<int> procs) const
{ return *this; }
2019-02-12 00:01:07 +05:00
};
2019-08-26 16:02:13 +05:00
2024-05-13 16:43:53 +05:00
inline void MyMPI_WaitAll (FlatArray<NG_MPI_Request> requests) { ; }
inline int MyMPI_WaitAny (FlatArray<NG_MPI_Request> requests) { return 0; }
2019-08-26 16:02:13 +05:00
2019-07-10 13:56:55 +05:00
#endif // PARALLEL
2019-02-12 00:01:07 +05:00
2019-08-27 19:52:54 +05:00
} // namespace ngcore
2019-02-12 00:01:07 +05:00
2019-07-10 13:56:55 +05:00
#endif // NGCORE_MPIWRAPPER_HPP
2019-02-12 00:01:07 +05:00