source: CIVL/examples/messagePassing/mpi.cvl@ 20a83c7

1.23 2.0 main test-branch
Last change on this file since 20a83c7 was 70e7b0f, checked in by Stephen Siegel <siegel@…>, 12 years ago

Separated out the point-to-point from the collective state by
creating a pair of comms for every MPI communicator.

git-svn-id: svn://vsl.cis.udel.edu/civl/trunk@778 fb995dde-84ed-4084-dfe6-e5aef3e2452c

  • Property mode set to 100644
File size: 3.9 KB
Line 
1
2typedef struct CMPI_Gcomm {
3 $gcomm p2p; // point-to-point communication
4 $gcomm col; // collective communication
5} CMPI_Gcomm;
6
7CMPI_Gcomm CMPI_Gcomm_create($scope scope, int size) {
8 CMPI_Gcomm result;
9
10 result.p2p = $gcomm_create(scope, size);
11 result.col = $gcomm_create(scope, size);
12 return result;
13}
14
15void CMPI_Gcomm_destroy(CMPI_Gcomm gc) {
16 $gcomm_destroy(gc.p2p);
17 $gcomm_destroy(gc.col);
18}
19
20typedef struct MPI_Comm {
21 $comm p2p; // point-to-point communication
22 $comm col; // collective communication
23} MPI_Comm;
24
25MPI_Comm MPI_Comm_create($scope scope, CMPI_Gcomm gc, int rank) {
26 MPI_Comm result;
27
28 result.p2p = $comm_create(scope, gc.p2p, rank);
29 result.col = $comm_create(scope, gc.col, rank);
30 return result;
31}
32
33void MPI_Comm_destroy(MPI_Comm comm) {
34 $comm_destroy(comm.p2p);
35 $comm_destroy(comm.col);
36}
37
38typedef enum {
39 MPI_INT,
40 MPI_FLOAT,
41 MPI_DOUBLE,
42 MPI_CHAR
43} MPI_Datatype;
44
45#define BCAST_TAG 999
46
47#define MPI_PROC_NULL -3
48
49typedef struct {
50 int MPI_SOURCE;
51 int MPI_TAG;
52 int MPI_ERROR;
53 int size;
54} MPI_Status;
55
56#define MPI_STATUS_IGNORE NULL
57
58#define MPI_STATUSES_IGNORE NULL
59
60int sizeofDatatype(MPI_Datatype datatype) {
61 switch (datatype) {
62 case MPI_INT:
63 return sizeof(int);
64 case MPI_FLOAT:
65 return sizeof(float);
66 case MPI_DOUBLE:
67 return sizeof(double);
68 case MPI_CHAR:
69 return sizeof(char);
70 default:
71 $assert(0, "Unreachable");
72 }
73}
74
75int MPI_Init() {
76 return 0;
77}
78
79int MPI_Finalize() {
80 return 0;
81}
82
83int MPI_Comm_size(MPI_Comm comm, int *size) {
84 *size = $comm_size(comm.p2p);
85 return 0;
86}
87
88int MPI_Comm_rank(MPI_Comm comm, int *rank) {
89 *rank = $comm_place(comm.p2p);
90 return 0;
91}
92
93
94int CMPI_Send(void *buf, int count, MPI_Datatype datatype, int dest,
95 int tag, $comm comm) {
96 if (dest >= 0) {
97 int size = count*sizeofDatatype(datatype);
98 int place = $comm_place(comm);
99 $message out = $message_pack(place, dest, tag, buf, size);
100
101 $comm_enqueue(comm, out);
102 }
103 return 0;
104}
105
106int MPI_Send(void *buf, int count, MPI_Datatype datatype, int dest,
107 int tag, MPI_Comm comm) {
108 return CMPI_Send(buf, count, datatype, dest, tag, comm.p2p);
109}
110
111
112int CMPI_Recv(void *buf, int count, MPI_Datatype datatype, int source,
113 int tag, $comm comm, MPI_Status *status) {
114 if (source >= 0) {
115 $message in = $comm_dequeue(comm, source, tag);
116 int size = count*sizeofDatatype(datatype);
117
118 $message_unpack(in, buf, size);
119 if (status != MPI_STATUS_IGNORE) {
120 status->size = $message_size(in);
121 status->MPI_SOURCE = $message_source(in);
122 status->MPI_TAG = $message_tag(in);
123 status->MPI_ERROR = 0;
124 }
125 }
126 return 0;
127}
128
129int MPI_Recv(void *buf, int count, MPI_Datatype datatype, int source,
130 int tag, MPI_Comm comm, MPI_Status *status) {
131 return CMPI_Recv(buf, count, datatype, source, tag, comm.p2p, status);
132}
133
134int MPI_Get_count(MPI_Status *status, MPI_Datatype datatype, int *count) {
135 *count = status->size/sizeofDatatype(datatype);
136 return 0;
137}
138
139int MPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
140 int dest, int sendtag,
141 void *recvbuf, int recvcount, MPI_Datatype recvtype,
142 int source, int recvtag,
143 MPI_Comm comm, MPI_Status *status) {
144 MPI_Send(sendbuf, sendcount, sendtype, dest, sendtag, comm);
145 MPI_Recv(recvbuf, recvcount, recvtype, source, recvtag, comm, status);
146 return 0;
147}
148
149/* Broadcasts a message from root to everyone else.
150 * Need to use a differnt comm.
151 */
152int MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
153 MPI_Comm comm) {
154 int place = $comm_place(comm.col);
155
156#ifdef DEBUG
157 printf("MPI_Bcast: place=%d, count=%d, root=%d\n", place, count, root);
158#endif
159 if (place == root) {
160 int nprocs = $comm_size(comm.col);
161
162 for (int i=0; i<nprocs; i++)
163 if (i != root)
164 CMPI_Send(buf, count, datatype, i, BCAST_TAG, comm.col);
165 } else {
166 CMPI_Recv(buf, count, datatype, root, BCAST_TAG, comm.col,
167 MPI_STATUS_IGNORE);
168 }
169 return 0;
170}
Note: See TracBrowser for help on using the repository browser.