source: CIVL/text/include/mpi.cvl@ d87ec9c

1.23 2.0 main test-branch
Last change on this file since d87ec9c was 42561ab, checked in by Stephen Siegel <siegel@…>, 12 years ago

Small changes to pack,unpack,operator application in meeting with Ziqing.

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

  • Property mode set to 100644
File size: 7.2 KB
Line 
1#ifndef __MPI_CVL__
2#define __MPI_CVL__
3//TODO make a Datatype struct, which has a field "int size;" Define one of these objects for MPI_INT, MPI_DOUBLE, etc.
4//TODO Then provide methods like MPI provides for creating new ones.
5//TODO then support MPI_Type_contig(datatype, int n).
6
7#define BCAST_TAG 999
8#define REDUCE_TAG 998
9
10/* Completed definition for mpi-common.h */
11struct MPI_Status{
12 int MPI_SOURCE;
13 int MPI_TAG;
14 int MPI_ERROR;
15 int size;
16};
17
18/* Definition of CIVL-MPI */
19typedef enum {
20 __UNINIT,
21 __INIT,
22 __FINALIZED
23}__MPI_Sys_status__;
24
25struct MPI_Request{
26 int id;
27};
28
29/* Definition of CMPI_Gcomm and MPI_Comm */
30typedef struct CMPI_Gcomm {
31 $gcomm p2p; // point-to-point communication
32 $gcomm col; // collective communication
33 $gbarrier gbarrier;
34} CMPI_Gcomm;
35
36struct MPI_Comm {
37 $comm p2p; // point-to-point communication
38 $comm col; // collective communication
39 $barrier barrier;
40 __MPI_Sys_status__ status;
41};
42
43/********************************** State **************************************/
44
45/* The number of times the MPI_Wtime function has been called */
46int CMPI_time_count = 0;
47
48/****************************** Helper Functions **********************************/
49int sizeofDatatype(MPI_Datatype datatype) {
50 switch (datatype) {
51 case MPI_INT:
52 return sizeof(int);
53 case MPI_FLOAT:
54 return sizeof(float);
55 case MPI_DOUBLE:
56 return sizeof(double);
57 case MPI_CHAR:
58 return sizeof(char);
59 default:
60 $assert(0, "Unreachable");
61 }
62}
63
64/************************** MPI LIB Implementations *******************************/
65
66$abstract double CMPI_time(int count);
67
68double MPI_Wtime() {
69 double result = CMPI_time(CMPI_time_count);
70
71 CMPI_time_count++;
72 return result;
73}
74
75CMPI_Gcomm CMPI_Gcomm_create($scope scope, int size) {
76 CMPI_Gcomm result;
77
78 result.p2p = $gcomm_create(scope, size);
79 result.col = $gcomm_create(scope, size);
80 result.gbarrier = $gbarrier_create(scope, size);
81 return result;
82}
83
84void CMPI_Gcomm_destroy(CMPI_Gcomm gc) {
85 $gcomm_destroy(gc.p2p);
86 $gcomm_destroy(gc.col);
87 $gbarrier_destroy(gc.gbarrier);
88}
89
90MPI_Comm CMPI_Comm_create($scope scope, CMPI_Gcomm gc, int rank) {
91 MPI_Comm result;
92
93 result.p2p = $comm_create(scope, gc.p2p, rank);
94 result.col = $comm_create(scope, gc.col, rank);
95 result.barrier = $barrier_create(scope, gc.gbarrier, rank);
96 result.status = __UNINIT;
97 return result;
98}
99
100void CMPI_Comm_destroy(MPI_Comm comm) {
101 $comm_destroy(comm.p2p);
102 $comm_destroy(comm.col);
103 $barrier_destroy(comm.barrier);
104}
105
106int __MPI_Init(MPI_Comm *comm) {
107 comm->status = __INIT;
108 return 0;
109}
110
111int __MPI_Finalize(MPI_Comm *comm) {
112 comm->status = __FINALIZED;
113 return 0;
114}
115
116int MPI_Comm_size(MPI_Comm comm, int *size) {
117 $assert(comm.status == __INIT, "MPI_Comm_size() cannot be invoked without MPI_Init() being called before.\n");
118 *size = $comm_size(comm.p2p);
119 return 0;
120}
121
122int MPI_Comm_rank(MPI_Comm comm, int *rank) {
123 $assert(comm.status == __INIT, "MPI_Comm_rank() cannot be invoked without MPI_Init() being called before.\n");
124 *rank = $comm_place(comm.p2p);
125 return 0;
126}
127
128
129int CMPI_Send(void *buf, int count, MPI_Datatype datatype, int dest,
130 int tag, $comm comm) {
131 if (dest >= 0) {
132 int size = count*sizeofDatatype(datatype);
133 int place = $comm_place(comm);
134 $message out = $message_pack(place, dest, tag, buf, size);
135
136 $comm_enqueue(comm, out);
137 }
138 return 0;
139}
140
141int MPI_Send(void *buf, int count, MPI_Datatype datatype, int dest,
142 int tag, MPI_Comm comm) {
143 $assert(comm.status == __INIT, "MPI_Send() cannot be invoked without MPI_Init() being called before.\n");
144 return CMPI_Send(buf, count, datatype, dest, tag, comm.p2p);
145}
146
147
148int CMPI_Recv(void *buf, int count, MPI_Datatype datatype, int source,
149 int tag, $comm comm, MPI_Status *status) {
150 if (source >= 0 || source == MPI_ANY_SOURCE) {
151 $message in = $comm_dequeue(comm, source, tag);
152 int size = count*sizeofDatatype(datatype);
153
154 $message_unpack(in, buf, size);
155 if (status != MPI_STATUS_IGNORE) {
156 status->size = $message_size(in);
157 status->MPI_SOURCE = $message_source(in);
158 status->MPI_TAG = $message_tag(in);
159 status->MPI_ERROR = 0;
160 }
161 }
162 return 0;
163}
164
165int MPI_Recv(void *buf, int count, MPI_Datatype datatype, int source,
166 int tag, MPI_Comm comm, MPI_Status *status) {
167 $assert(comm.status == __INIT,
168 "MPI_Recv() cannot be invoked without "
169 "MPI_Init() being called before.\n");
170 return CMPI_Recv(buf, count, datatype, source, tag, comm.p2p, status);
171}
172
173int MPI_Get_count(MPI_Status *status, MPI_Datatype datatype, int *count) {
174 *count = status->size/sizeofDatatype(datatype);
175 return 0;
176}
177
178int MPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
179 int dest, int sendtag,
180 void *recvbuf, int recvcount, MPI_Datatype recvtype,
181 int source, int recvtag,
182 MPI_Comm comm, MPI_Status *status) {
183 $assert(comm.status == __INIT,
184 "MPI_Sendrecv() cannot be invoked "
185 "without MPI_Init() being called before.\n");
186 // not correct for checking potential deadlock...rewrite:
187 MPI_Send(sendbuf, sendcount, sendtype, dest, sendtag, comm);
188 MPI_Recv(recvbuf, recvcount, recvtype, source, recvtag, comm, status);
189 return 0;
190}
191
192/* Broadcasts a message from root to everyone else.
193 * Need to use a differnt comm.
194 */
195int MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
196 MPI_Comm comm) {
197 $assert(comm.status == __INIT,
198 "MPI_Bcast() cannot be invoked without MPI_Init() "
199 "being called before.\n");
200 if ($comm_place(comm.col) == root) {
201 int nprocs = $comm_size(comm.col);
202
203 for (int i=0; i<nprocs; i++)
204 if (i != root)
205 CMPI_Send(buf, count, datatype, i, BCAST_TAG, comm.col);
206 } else {
207 CMPI_Recv(buf, count, datatype, root, BCAST_TAG, comm.col,
208 MPI_STATUS_IGNORE);
209 }
210 return 0;
211}
212
213/* Reduces values on all processes to a single value */
214int MPI_Reduce(void* sendbuf, void* recvbuf, int count,
215 MPI_Datatype datatype, MPI_Op op, int root,
216 MPI_Comm comm) {
217 $assert(comm.status == __INIT,
218 "MPI_Reduce() cannot be invoked without "
219 "MPI_Init() being called before.\n");
220 CMPI_Send(sendbuf, count, datatype, root, REDUCE_TAG, comm.col);
221 if ($comm_place(comm.col) == root) {
222 int nprocs = $comm_size(comm.col);
223
224 for (int i = 0; i<nprocs; i++) {
225 $message in = $comm_dequeue(comm.col, i, REDUCE_TAG);
226 int size = count * sizeofDatatype(datatype);
227
228 $bundle_unpack_apply(in.data, recvbuf, count, op);
229 $assert(in.size <= size,
230 "Message of size %d exceeds the specified size %d.",
231 in.size, size);
232 }
233 }
234 return 0;
235}
236
237/* Combines values from all processes and distributes the result back to all processes */
238/* default root is 0 */
239int MPI_Allreduce(void* sendbuf, void* recvbuf, int count,
240 MPI_Datatype datatype,
241 MPI_Op op, MPI_Comm comm) {
242 int root = 0;
243 MPI_Status status;
244
245 $assert(comm.status == __INIT,
246 "MPI_Allreduce() cannot be invoked without "
247 "MPI_Init() being called before.\n");
248 MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm);
249 MPI_Bcast(recvbuf, count, datatype, root, comm);
250 return 0;
251}
252
253int MPI_Barrier(MPI_Comm comm){
254
255 $assert(comm.status == __INIT, "MPI_Allreduce() cannot be invoked without MPI_Init() being called before.\n");
256 $barrier_call(comm.barrier);
257}
258#endif
Note: See TracBrowser for help on using the repository browser.