Fix warnings and errors in ampi_noimpl.C
[charm.git] / src / util / pup_mpi.h
blob9d78ff6f84d204f0b73be3c9b81c64e705520dfd
1 /*
2 PUP -> MPI interface routines
5 Orion Sky Lawlor, olawlor@acm.org, 2004/9/15
6 */
7 #ifndef __UIUC_CHARM_PUP_MPI_H
8 #define __UIUC_CHARM_PUP_MPI_H
10 #include "pup.h"
11 #include "mpi.h"
13 #define pup_checkMPI(err) pup_checkMPIerr(err,__FILE__,__LINE__);
14 inline void pup_checkMPIerr(int mpi_err,const char *file,int line) {
15 if (mpi_err!=MPI_SUCCESS) {
16 CmiError("MPI Routine returned error %d at %s:%d\n",
17 mpi_err,file,line);
18 CmiAbort("MPI Routine returned error code");
22 /// Return the number of dt's in the next MPI message from/tag/comm.
23 inline int MPI_Incoming_pup(MPI_Datatype dt,int from,int tag,MPI_Comm comm) {
24 MPI_Status sts;
25 pup_checkMPI(MPI_Probe(from,tag,comm,&sts));
26 int len; pup_checkMPI(MPI_Get_count(&sts,dt,&len));
27 return len;
30 /// MPI_Recv, but using an object T with a pup routine.
31 template <class T>
32 inline void MPI_Recv_pup(T &t, int from,int tag,MPI_Comm comm) {
33 int len=MPI_Incoming_pup(MPI_BYTE,from,tag,comm);
34 MPI_Status sts;
35 char *buf=new char[len];
36 pup_checkMPI(MPI_Recv(buf,len,MPI_BYTE, from,tag,comm,&sts));
37 PUP::fromMemBuf(t,buf,len);
38 delete[] buf;
41 /// MPI_Send, but using an object T with a pup routine.
42 template <class T>
43 inline void MPI_Send_pup(T &t, int to,int tag,MPI_Comm comm) {
44 size_t len=PUP::size(t); char *buf=new char[len];
45 PUP::toMemBuf(t,buf,len);
46 pup_checkMPI(MPI_Send(buf,len,MPI_BYTE, to,tag,comm));
47 delete[] buf;
50 /// MPI_Bcast, but using an object T with a pup routine.
51 template <class T>
52 inline void MPI_Bcast_pup(T &t, int root,MPI_Comm comm) {
53 int myRank;
54 MPI_Comm_rank(comm,&myRank);
55 /* Can't do broadcast until everybody knows the size */
56 size_t len=0;
57 if(myRank == root) len=PUP::size(t);
58 pup_checkMPI(MPI_Bcast(&len,1,MPI_INT,root,comm));
59 /* Now pack object and send off */
60 char *buf=new char[len];
61 if(myRank == root) PUP::toMemBuf(t,buf,len);
62 pup_checkMPI(MPI_Bcast(buf,len,MPI_BYTE, root,comm));
63 PUP::fromMemBuf(t,buf,len);
64 delete [] buf;
68 #endif