source: CIVL/include/impls/mpi.cvl@ 1aaefd4

main test-branch
Last change on this file since 1aaefd4 was ea777aa, checked in by Alex Wilton <awilton@…>, 3 years ago

Moved examples, include, build_default.properties, common.xml, and README out from dev.civl.com into the root of the repo.

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

  • Property mode set to 100644
File size: 25.1 KB
Line 
1#ifndef __CIVL_MPI__
2#define __CIVL_MPI__
3
4#include <seq.cvh>
5#include <civlc.cvh>
6#include <concurrency.cvh>
7#include <bundle.cvh>
8#include <mpi.h>
9#include <civl-mpi.cvh>
10#include <stdlib.h>
11#include <string.h>
12#include <pointer.cvh>
13#include <collate.cvh>
14
15
16extern const int MPI_IN_PLACE_SPOT = 0;
17
18/* Completed definition for mpi-common.h */
19
20$mpi_state _mpi_state=_MPI_UNINIT;
21
22/************************** MPI LIB Implementations *******************************/
23
24int $mpi_init(void) {
25 $assert(_mpi_state == _MPI_UNINIT, "Process can only call MPI_Init() at most once.");
26 _mpi_state = _MPI_INIT;
27 return 0;
28}
29
30int MPI_Finalize(void) {
31 $assert(_mpi_state == _MPI_INIT, "Process can only call MPI_Finalize() after the "
32 "MPI enviroment is created and before cleaned.");
33 _mpi_state = _MPI_FINALIZED;
34 return 0;
35}
36
37double MPI_Wtime() {
38 double result;
39 int CMPI_time_count = $next_time_count();
40
41 $assert(_mpi_state == _MPI_INIT, "MPI_Wtime() cannot be invoked "
42 "without MPI_Init() being called before.\n");
43 result = $mpi_time(CMPI_time_count);
44 if (CMPI_time_count > 0) {
45 $assume(result > $mpi_time(CMPI_time_count-1));
46 } else {
47 $assume(result > 0);
48 }
49 return result;
50}
51
52int MPI_Comm_size(MPI_Comm comm, int *size) {
53#ifndef _MPI_CONTRACT
54 $assert(_mpi_state == _MPI_INIT, "MPI_Comm_size() cannot be "
55 "invoked without MPI_Init() being called before.\n");
56#endif
57 *size = $mpi_comm_size(comm);
58 return 0;
59}
60
61int MPI_Comm_rank(MPI_Comm comm, int *rank) {
62#ifndef _MPI_CONTRACT
63 $assert(_mpi_state == _MPI_INIT, "MPI_Comm_rank() cannot be "
64 "invoked without MPI_Init() being called before.\n");
65#endif
66 *rank = $mpi_comm_place(comm);
67 return 0;
68}
69
70int MPI_Send(const void *buf, int count, MPI_Datatype datatype, int dest,
71 int tag, MPI_Comm comm) {
72 $assert(_mpi_state == _MPI_INIT, "MPI_Send() cannot be invoked "
73 "without MPI_Init() being called before.\n");
74
75#ifdef _MPI_NON_BLOCKING
76 MPI_Request request;
77
78 $mpi_isend(buf, count, datatype, dest, tag, comm, &request);
79 $mpi_wait(&request, MPI_STATUS_IGNORE);
80#elif defined(_MPI_BLOCKING)
81 $mpi_send(buf, count, datatype, dest, tag, comm);
82#endif
83 return 1;
84}
85
86int MPI_Isend(void *buf, int count, MPI_Datatype datatype, int dest,
87 int tag, MPI_Comm comm, MPI_Request * request) {
88 $assert(_mpi_state == _MPI_INIT, "MPI_Isend() cannot be invoked "
89 "without MPI_Init() being called before.\n");
90#ifdef _MPI_NON_BLOCKING
91 return $mpi_isend(buf, count, datatype, dest, tag, comm, request);
92#elif defined(_MPI_BLOCKING)
93 $assert(0, "MPI_Isend is not supported in the MPI blocking implementation.");
94 return 0;
95#endif
96}
97
98int MPI_Recv(void *buf, int count, MPI_Datatype datatype, int source,
99 int tag, MPI_Comm comm, MPI_Status *status) {
100 $assert(_mpi_state == _MPI_INIT, "MPI_Recv() cannot be invoked "
101 "without MPI_Init() being called before.\n");
102#ifdef _MPI_NON_BLOCKING
103 MPI_Request request;
104
105 $mpi_irecv(buf, count, datatype, source, tag, comm, &request);
106 $mpi_wait(&request, status);
107#elif defined(_MPI_BLOCKING)
108 $mpi_recv(buf, count, datatype, source, tag, comm, status);
109#endif
110 return 1;
111}
112
113int MPI_Irecv(void *buf, int count, MPI_Datatype datatype, int source,
114 int tag, MPI_Comm comm, MPI_Request * request) {
115 $assert(_mpi_state == _MPI_INIT, "MPI_Irecv() cannot be invoked "
116 "without MPI_Init() being called before.\n");
117#ifdef _MPI_NON_BLOCKING
118 return $mpi_irecv(buf, count, datatype, source, tag, comm, request);
119#elif defined(_MPI_BLOCKING)
120 $assert(0, "MPI_Irecv is not supported in the MPI blocking implementation.");
121 return 0;
122#endif
123}
124
125int MPI_Wait(MPI_Request * request, MPI_Status * status) {
126 $assert(_mpi_state == _MPI_INIT, "MPI_Wait() cannot be invoked "
127 "without MPI_Init() being called before.\n");
128#ifdef _MPI_NON_BLOCKING
129 $mpi_wait(request, status);
130 return 1;
131#elif defined(_MPI_BLOCKING)
132 $assert(0, "MPI_Wait is not supported in the MPI blocking implementation.");
133 return 0;
134#endif
135}
136
137int MPI_Waitall(int count, MPI_Request array_of_requests[], MPI_Status array_of_statuses[]) {
138 $assert(_mpi_state == _MPI_INIT, "MPI_Waitall() cannot be invoked "
139 "without MPI_Init() being called before.\n");
140#ifdef _MPI_NON_BLOCKING
141 $for (int i : 0 .. count-1) {
142 MPI_Status * status = array_of_statuses == MPI_STATUSES_IGNORE ? MPI_STATUS_IGNORE : array_of_statuses + i;
143 MPI_Request * req = array_of_requests + i;
144
145 $mpi_wait(req, status);
146 }
147 return 1;
148#elif defined(_MPI_BLOCKING)
149 $assert(0, "MPI_Waitall is not supported in the MPI blocking implementation.");
150 return 0;
151#endif
152}
153
154int MPI_Test(MPI_Request *request, int *flag, MPI_Status *status) {
155 $assert(_mpi_state == _MPI_INIT, "MPI_Test() cannot be invoked "
156 "without MPI_Init() being called before.\n");
157#ifdef _MPI_NON_BLOCKING
158 if (*request == MPI_REQUEST_NULL) {
159 *flag = 1;
160 return 1;
161 }
162 $choose {
163 $when ($true) {
164 $mpi_wait(request, status);
165 *flag = 1;
166 }
167 $when ($true)
168 *flag = 0;
169 }
170 return 1;
171#elif defined(_MPI_BLOCKING)
172 $assert(0, "MPI_Test is not supported in the MPI blocking implementation.");
173 return 0;
174#endif
175}
176
177int MPI_Request_free(MPI_Request * request) {
178 // The standard does not say this function accepts MPI_REQUEST_NULL
179#ifdef _MPI_NON_BLOCKING
180 $free(*request);
181 *request = MPI_REQUEST_NULL;
182 return 1;
183#elif defined(_MPI_BLOCKING)
184 $assert(0, "MPI_Request_free is not supported in the MPI blocking implementation.");
185 return 0;
186#endif
187}
188
189int MPI_Get_count(MPI_Status *status, MPI_Datatype datatype, int *count) {
190#ifndef _MPI_CONTRACT
191 $assert(_mpi_state == _MPI_INIT, "MPI_Get_count() cannot be invoked "
192 "without MPI_Init() being called before.\n");
193#endif
194 *count = status->size/sizeofDatatype(datatype);
195 return 0;
196}
197
198int MPI_Get_processor_name(char * name, int * resultlen) {
199 $abstract int MPI_GET_PROCESSOR_NAME(char *, int *);
200
201 return MPI_GET_PROCESSOR_NAME(name, resultlen);
202}
203
204int MPI_Sendrecv(void *sendbuf, int sendcount, MPI_Datatype sendtype,
205 int dest, int sendtag,
206 void *recvbuf, int recvcount, MPI_Datatype recvtype,
207 int source, int recvtag,
208 MPI_Comm comm, MPI_Status *status) {
209 $assert(_mpi_state == _MPI_INIT, "MPI_Sendrecv() cannot be invoked "
210 "without MPI_Init() being called before.\n");
211#ifdef _MPI_CONTRACT
212 $elaborate(dest);
213 $elaborate(source);
214#else
215 $mpi_check_buffer(sendbuf, sendcount, sendtype);
216#endif
217 // not correct for checking potential deadlock...rewrite:
218 $mpi_sendrecv(sendbuf, sendcount, sendtype, dest, sendtag, recvbuf, recvcount,
219 recvtype, source, recvtag, comm, status);
220 return 0;
221}
222/******************************** Collective ***********************************/
223/* Broadcasts a message from root to everyone else.
224 * Need to use a differnt comm.
225 */
226int MPI_Bcast(void *buf, int count, MPI_Datatype datatype, int root,
227 MPI_Comm comm) {
228 int place = $mpi_comm_place(comm);
229 int nprocs = $mpi_comm_size(comm);
230 int datatypes[1] = {(int)datatype};
231 // MPI library defined collective operation checking entries:
232 $bundle checkerEntry; //the checking entry of this call
233 $bundle specEntry; //a recorded entry as specification
234
235#ifndef _MPI_CONTRACT
236 $assert (_mpi_state == _MPI_INIT,
237 "MPI_Bcast() cannot be invoked without MPI_Init() "
238 "being called before.\n");
239#endif
240 if(place == root)
241 $mpi_check_buffer(buf, count, datatype);
242 checkerEntry = $mpi_create_coroutine_entry(BCAST_TAG, root, -1, 1, datatypes);
243 specEntry = $mpi_check_collective(comm, checkerEntry);
244 $mpi_diff_coroutine_entries(specEntry, checkerEntry, place);
245 $mpi_bcast(buf, count, datatype, root, BCAST_TAG, comm, "MPI_Bcast()");
246 return 0;
247}
248
249/* Reduces values on all processes to a single value */
250int MPI_Reduce(const void* sendbuf, void* recvbuf, int count,
251 MPI_Datatype datatype, MPI_Op op, int root,
252 MPI_Comm comm) {
253 int place = $mpi_comm_place(comm);
254 int nprocs = $mpi_comm_size(comm);
255 int datatypes[1] = {(int)datatype};
256 // MPI library defined collective operation checking entries:
257 $bundle checkerEntry; //the checking entry of this call
258 $bundle specEntry; //a recorded entry as specification
259
260#ifndef _MPI_CONTRACT
261 $assert (_mpi_state == _MPI_INIT,
262 "MPI_Reduce() cannot be invoked without "
263 "MPI_Init() being called before.\n");
264#endif
265 checkerEntry = $mpi_create_coroutine_entry(REDUCE_TAG, root, (int)op, 1,
266 datatypes);
267 specEntry = $mpi_check_collective(comm, checkerEntry);
268 $mpi_diff_coroutine_entries(specEntry, checkerEntry, place);
269 $mpi_check_buffer(sendbuf, count, datatype);
270 $assert(0 <= op && op <= 13, "unknown MPI reduce operation"); // refer to op.h & mpi.h in ABC/src/include for how MPI_Op is defined
271 $mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, REDUCE_TAG, comm, "MPI_Reduce()");
272 return 0;
273}
274
275/* Combines values from all processes and distributes the result back to all processes */
276/* default root is 0 */
277int MPI_Allreduce(const void* sendbuf, void* recvbuf, int count,
278 MPI_Datatype datatype,
279 MPI_Op op, MPI_Comm comm) {
280 int root = 0;
281 int place = $mpi_comm_place(comm);
282 int nprocs = $mpi_comm_size(comm);
283 int datatypes[1] = {(int)datatype};
284 MPI_Status status;
285 // MPI library defined collective operation checking entries:
286 $bundle checkerEntry; //the checking entry of this call
287 $bundle specEntry; //a recorded entry as specification
288
289#ifndef _MPI_CONTRACT
290 $assert(_mpi_state == _MPI_INIT, "MPI_Allreduce() cannot be invoked without "
291 "MPI_Init() being called before.\n");
292#endif
293 $mpi_check_buffer(sendbuf, count, datatype);
294 checkerEntry = $mpi_create_coroutine_entry(ALLREDUCE_TAG, root,
295 (int)op, 1, datatypes);
296 specEntry = $mpi_check_collective(comm, checkerEntry);
297 $mpi_diff_coroutine_entries(specEntry, checkerEntry, place);
298 $assert(0 <= op && op <= 13, "unknown MPI reduce operation"); // refer to op.h & mpi.h in ABC/src/include for how MPI_Op is defined
299 $mpi_reduce(sendbuf, recvbuf, count, datatype, op, root, ALLREDUCE_TAG, comm,
300 "MPI_Allreduce()");
301 $mpi_bcast(recvbuf, count, datatype, root, ALLREDUCE_TAG, comm,
302 "MPI_Allreduce()");
303 return 0;
304}
305
306int MPI_Barrier(MPI_Comm comm){
307 int place = $mpi_comm_place(comm);
308 int nprocs = $mpi_comm_size(comm);
309 // MPI library defined collective operation checking entries:
310 $bundle checkerEntry; //the checking entry of this call
311 $bundle specEntry; //a recorded entry as specification
312
313#ifndef _MPI_CONTRACT
314 $assert(_mpi_state == _MPI_INIT, "MPI_Barrier() cannot be invoked "
315 "without MPI_Init() being called before.\n");
316#endif
317 checkerEntry = $mpi_create_coroutine_entry(BARRIER_TAG, 0, -1,
318 0, NULL);
319 specEntry = $mpi_check_collective(comm, checkerEntry);
320 $mpi_diff_coroutine_entries(specEntry, checkerEntry, place);
321 $mpi_barrier(comm);
322 return 0;
323}
324
325/* 1. If comm is an intracommunicator, each process (includes root process) sends the content
326 of its send buffer to the root process. Root process receives the messages and stores
327 them in rank order
328 2. TODO: If comm is an intercommunicator, it's not supported yet */
329int MPI_Gather(const void* sendbuf, int sendcount, MPI_Datatype sendtype,
330 void* recvbuf, int recvcount, MPI_Datatype recvtype,
331 int root, MPI_Comm comm){
332 int place = $mpi_comm_place(comm);
333 int nprocs = $mpi_comm_size(comm);
334 int datatypes[2] = {(int)sendtype, (int)recvtype};
335 // MPI library defined collective operation checking entries:
336 $bundle checkerEntry; //the checking entry of this call
337 $bundle specEntry; //a recorded entry as specification
338
339#ifndef _MPI_CONTRACT
340 $assert(_mpi_state == _MPI_INIT, "MPI_Gather() cannot be invoked without "
341 "MPI_Init() being called before.\n");
342#endif
343 if(sendbuf != MPI_IN_PLACE)
344 $mpi_check_buffer(sendbuf, sendcount, sendtype);
345 checkerEntry = $mpi_create_coroutine_entry(GATHER_TAG, root, -1, 2, datatypes);
346 specEntry = $mpi_check_collective(comm, checkerEntry);
347 $mpi_diff_coroutine_entries(specEntry, checkerEntry, place);
348 $mpi_gather(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype,
349 root, GATHER_TAG, comm, "MPI_Gather()");
350 return 0;
351}
352
353/* The inverse operation of MPI_Gather() */
354int MPI_Scatter(const void* sendbuf, int sendcount, MPI_Datatype sendtype,
355 void* recvbuf, int recvcount, MPI_Datatype recvtype, int root,
356 MPI_Comm comm){
357 int place = $mpi_comm_place(comm);
358 int nprocs = $mpi_comm_size(comm);
359 int datatypes[2] = {(int)sendtype, (int)recvtype};
360 // MPI library defined collective operation checking entries:
361 $bundle checkerEntry; //the checking entry of this call
362 $bundle specEntry; //a recorded entry as specification
363
364#ifndef _MPI_CONTRACT
365 $assert(_mpi_state == _MPI_INIT, "MPI_Scatter() cannot be invoked without "
366 "MPI_Init() being called before.\n");
367#endif
368 if (place == root)
369 $mpi_check_buffer(sendbuf, sendcount, sendtype);
370 checkerEntry = $mpi_create_coroutine_entry(SCATTER_TAG, root, -1, 2,
371 datatypes);
372 specEntry = $mpi_check_collective(comm, checkerEntry);
373 $mpi_diff_coroutine_entries(specEntry, checkerEntry, place);
374 $mpi_scatter(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype,
375 root, SCATTER_TAG, comm, "MPI_Scatter()");
376 return 0;
377}
378
379
380/* MPI_Gatherv extends the functionality of MPI_Gather by allowing a varying count of data to be sent to root process, since recvcounts is now an array.*/
381int MPI_Gatherv(const void* sendbuf, int sendcount, MPI_Datatype sendtype,
382 void* recvbuf, const int recvcounts[], const int displs[],
383 MPI_Datatype recvtype, int root, MPI_Comm comm){
384 int place = $mpi_comm_place(comm);
385 int nprocs = $mpi_comm_size(comm);
386 int datatypes[2] = {(int)sendtype, (int)recvtype};
387 int recvcount = 0;
388 // MPI library defined collective operation checking entries:
389 $bundle checkerEntry; //the checking entry of this call
390 $bundle specEntry; //a recorded entry as specification
391
392#ifndef _MPI_CONTRACT
393 $assert(_mpi_state == _MPI_INIT, "MPI_Gatherv() cannot be invoked without "
394 "MPI_Init() being called before.\n");
395#endif
396 if(sendbuf != MPI_IN_PLACE)
397 $mpi_check_buffer(sendbuf, sendcount, sendtype);
398 checkerEntry = $mpi_create_coroutine_entry(GATHERV_TAG, root,
399 -1, 2, datatypes);
400 specEntry = $mpi_check_collective(comm, checkerEntry);
401 $mpi_diff_coroutine_entries(specEntry, checkerEntry, place);
402 $mpi_gatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs,
403 recvtype, root, GATHERV_TAG, comm, "MPI_Gatherv()");
404 return 0;
405}
406
407/* MPI_Scatterv extends the functionality of MPI_Scatter by allowing a varying count of data to be sent to each process, since sendcounts is now an array.*/
408int MPI_Scatterv(const void* sendbuf, const int sendcounts[], const
409 int displs[], MPI_Datatype sendtype, void* recvbuf,
410 int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm){
411 int place = $mpi_comm_place(comm);
412 int nprocs = $mpi_comm_size(comm);
413 int datatypes[2] = {(int)sendtype, (int)recvtype};
414 int sendcount = 0;
415 // MPI library defined collective operation checking entries:
416 $bundle checkerEntry; //the checking entry of this call
417 $bundle specEntry; //a recorded entry as specification
418
419#ifndef _MPI_CONTRACT
420 $assert(_mpi_state == _MPI_INIT, "MPI_Scatterv() cannot be invoked without "
421 "MPI_Init() being called before.\n");
422#endif
423 if (place == root) {
424 for (int i = 0; i < nprocs; i++) sendcount += sendcounts[i];
425 $mpi_check_buffer(sendbuf, sendcount, sendtype);
426 }
427 checkerEntry = $mpi_create_coroutine_entry(SCATTERV_TAG,
428 root, -1, 2,
429 datatypes);
430 specEntry = $mpi_check_collective(comm, checkerEntry);
431 $mpi_diff_coroutine_entries(specEntry, checkerEntry, place);
432 $mpi_scatterv(sendbuf, sendcounts, displs, sendtype, recvbuf,
433 recvcount, recvtype, root, SCATTERV_TAG, comm,
434 "MPI_Scatterv()");
435 return 0;
436}
437
438int MPI_Allgather(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
439 void *recvbuf, int recvcount, MPI_Datatype recvtype,
440 MPI_Comm comm){
441 int place = $mpi_comm_place(comm);
442 int nprocs = $mpi_comm_size(comm);
443 int datatypes[2] = {(int)sendtype, (int)recvtype};
444 // MPI library defined collective operation checking entries:
445 $bundle checkerEntry; //the checking entry of this call
446 $bundle specEntry; //a recorded entry as specification
447
448#ifndef _MPI_CONTRACT
449 $assert(_mpi_state == _MPI_INIT, "MPI_Allgather() cannot be invoked without "
450 "MPI_Init() being called before.\n");
451#endif
452 if(sendbuf != MPI_IN_PLACE)
453 $mpi_check_buffer(sendbuf, sendcount, sendtype);
454 checkerEntry = $mpi_create_coroutine_entry(ALLGATHER_TAG, 0, -1, 2,
455 datatypes);
456 specEntry = $mpi_check_collective(comm, checkerEntry);
457 $mpi_diff_coroutine_entries(specEntry, checkerEntry, place);
458
459 if (sendbuf != MPI_IN_PLACE)
460 $mpi_gather(sendbuf, sendcount, sendtype,
461 recvbuf, recvcount, recvtype,
462 0, ALLGATHER_TAG, comm, "MPI_Allgather()");
463 else {
464 void * in_buf = $mpi_malloc(recvcount, recvtype);
465
466 memcpy(in_buf, recvbuf + recvcount*place, sizeofDatatype(recvtype) * recvcount);
467 $mpi_gather(in_buf, recvcount, recvtype,
468 recvbuf, recvcount, recvtype,
469 0, ALLGATHER_TAG, comm, "MPI_Allgather()");
470 $free(in_buf);
471 }
472 $mpi_bcast(recvbuf, recvcount*nprocs, recvtype, 0, ALLGATHER_TAG, comm,
473 "MPI_Allgather()");
474 return 0;
475}
476
477int MPI_Reduce_scatter(const void *sendbuf, void *recvbuf, const int recvcount[],
478 MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) {
479 int total_count, i;
480 int nprocs = $mpi_comm_size(comm);
481 int rank = $mpi_comm_place(comm);
482 int root = 0;
483 int displs[nprocs];
484 int datatypes[1] = {(int)datatype};
485
486 // MPI library defined collective operation checking entries:
487 $bundle checkerEntry; //the checking entry of this call
488 $bundle specEntry; //a recorded entry as specification
489
490#ifndef _MPI_CONTRACT
491 $assert(_mpi_state == _MPI_INIT, "MPI_Reduce_scatter() cannot be invoked without "
492 "MPI_Init() being called before.\n");
493#endif
494 $mpi_check_buffer(sendbuf, recvcount[rank], datatype);
495 for(total_count = 0, i = 0; i<nprocs; i++) {
496 displs[i] = total_count;
497 total_count += recvcount[i];
498 }
499 checkerEntry = $mpi_create_coroutine_entry(RED_SCATTER_TAG, root, (int)op, 1,
500 datatypes);
501 specEntry = $mpi_check_collective(comm, checkerEntry);
502 $mpi_diff_coroutine_entries(specEntry, checkerEntry, rank);
503 /* Note: In MPI standard, the sendbuf and recvbuf shall not be the
504 * same, the implementation here is a lower layer helper function
505 * for MPI_Reduce routine, and the reason it plays a trick here is
506 * because allocating a memory space for a void pointer is not
507 * allowed in CIVL yet. */
508 void * temp = $mpi_malloc(total_count, datatype);
509
510 $assert(0 <= op && op <= 13, "unknown MPI reduce operation"); // refer to op.h & mpi.h in ABC/src/include for how MPI_Op is defined
511 $mpi_reduce(sendbuf, temp, total_count, datatype, op,
512 root, RED_SCATTER_TAG, comm, "MPI_Reduce_scatter()");
513 $mpi_scatterv(temp, recvcount, displs, datatype, recvbuf,
514 recvcount[rank], datatype, root, RED_SCATTER_TAG, comm,
515 "MPI_Reduce_scatter()");
516 free(temp);
517 return 0;
518}
519
520int MPI_Alltoall(const void *sendbuf, int sendcount, MPI_Datatype sendtype,
521 void *recvbuf, int recvcount, MPI_Datatype recvtype,
522 MPI_Comm comm) {
523 int nprocs = $mpi_comm_size(comm);
524 int rank = $mpi_comm_place(comm);
525 int root = 0;
526 int displs[nprocs];
527 int sendcounts[nprocs];
528 int datatypes[2] = {(int)sendtype, (int)recvtype};
529
530 // MPI library defined collective operation checking entries:
531 $bundle checkerEntry; //the checking entry of this call
532 $bundle specEntry; //a recorded entry as specification
533
534 for(int i=0; i<nprocs; i++) {
535 sendcounts[i] = sendcount;
536 displs[i] = (i == 0)? 0 : (displs[i-1] + sendcount);
537 }
538#ifndef _MPI_CONTRACT
539 $assert(_mpi_state == _MPI_INIT, "MPI_Alltoall() cannot be invoked without "
540 "MPI_Init() being called before.\n");
541#endif
542 $mpi_check_buffer(sendbuf, sendcount * nprocs, sendtype);
543 checkerEntry = $mpi_create_coroutine_entry(ALLTOALL_TAG, root, -1, 2,
544 datatypes);
545 specEntry = $mpi_check_collective(comm, checkerEntry);
546 $mpi_diff_coroutine_entries(specEntry, checkerEntry, rank);
547 for(int i = 0; i < nprocs; i++) {
548 void * ptr = $mpi_pointer_add(recvbuf, i*sendcount, recvtype);
549
550 $mpi_scatterv(sendbuf, sendcounts, displs, sendtype,
551 ptr, recvcount, recvtype, i, ALLTOALL_TAG, comm,
552 "MPI_Alltoall()");
553 }
554 return 0;
555}
556
557int MPI_Alltoallv(const void* sendbuf, const int sendcounts[],
558 const int sdispls[], MPI_Datatype sendtype, void* recvbuf,
559 const int recvcounts[], const int rdispls[], MPI_Datatype recvtype,
560 MPI_Comm comm) {
561 int nprocs = $mpi_comm_size(comm);
562 int place = $mpi_comm_place(comm);
563 int datatypes[2] = {(int)sendtype, (int)recvtype};
564 int sendcount = 0;
565 int recvcount = 0;
566 // MPI library defined collective operation checking entries:
567 $bundle checkerEntry; //the checking entry of this call
568 $bundle specEntry; //a recorded entry as specification
569
570#ifndef _MPI_CONTRACT
571 $assert(_mpi_state == _MPI_INIT, "MPI_Alltoallv() cannot be invoked without "
572 "MPI_Init() being called before.\n");
573#endif
574 for (int i = 0; i < nprocs; i++) {
575 sendcount += sendcounts[i];
576 recvcount += recvcounts[i];
577 }
578 $mpi_check_buffer(sendbuf, sendcount, sendtype);
579 checkerEntry = $mpi_create_coroutine_entry(ALLTOALLV_TAG, 0, -1, 2,
580 datatypes);
581 specEntry = $mpi_check_collective(comm, checkerEntry);
582 $mpi_diff_coroutine_entries(specEntry, checkerEntry, place);
583 for(int i = 0; i < nprocs; i++) {
584 void * ptr = $mpi_pointer_add(recvbuf, rdispls[i], recvtype);
585
586 $mpi_scatterv(sendbuf, sendcounts, sdispls, sendtype,
587 ptr, recvcounts[i], recvtype, i, ALLTOALLV_TAG, comm,
588 "MPI_Alltoallv()");
589 }
590 return 0;
591}
592
593int MPI_Alltoallw(const void* sendbuf, const int sendcounts[], const int sdispls[],
594 const MPI_Datatype sendtypes[], void* recvbuf,
595 const int recvcounts[], const int rdispls[],
596 const MPI_Datatype recvtypes[], MPI_Comm comm) {
597 int nprocs = $mpi_comm_size(comm);
598 int place = $mpi_comm_place(comm);
599 int sdispls_offset[nprocs];
600
601 for (int i = 0; i < nprocs; i++)
602 sdispls_offset[i] = sdispls[i] / sizeofDatatype(sendtypes[i]);
603#ifndef _MPI_CONTRACT
604 $assert(_mpi_state == _MPI_INIT, "MPI_Alltoallw() cannot be invoked without "
605 "MPI_Init() being called before.\n");
606#endif
607 for(int i = 0; i < nprocs; i++) {
608 int recv_t_size = sizeofDatatype(recvtypes[i]);
609 void * ptr = $mpi_pointer_add(recvbuf, rdispls[i] / recv_t_size, recvtypes[i]);
610 void * sendptr = $mpi_pointer_add(sendbuf, sdispls_offset[i], sendtypes[i]);
611
612 $mpi_check_buffer(sendptr, sendcounts[i], sendtypes[i]);
613 $mpi_scatterv(sendbuf, sendcounts, sdispls_offset, sendtypes[i],
614 ptr, recvcounts[i], recvtypes[i], i,
615 ALLTOALLW_TAG, comm, "MPI_Alltoallw()");
616 }
617 return 0;
618}
619
620int MPI_Scan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
621 MPI_Op op, MPI_Comm comm) {
622#ifndef _MPI_CONTRACT
623 $assert(_mpi_state == _MPI_INIT, "MPI_Scan() cannot be invoked without "
624 "MPI_Init() being called before.\n");
625#endif
626 int place = $mpi_comm_place(comm);
627 int datatype_enum2int = (int)datatype;
628
629 // check consistency of a group of collective routine calls
630 $bundle checkerEntry = $mpi_create_coroutine_entry(SCAN_TAG, -1, -1, 1,
631 &datatype_enum2int);
632 $bundle specEntry = $mpi_check_collective(comm, checkerEntry);
633
634 $mpi_diff_coroutine_entries(specEntry, checkerEntry, place);
635 if (sendbuf != MPI_IN_PLACE)
636 $mpi_check_buffer(sendbuf, count, datatype);
637 $mpi_check_buffer(recvbuf, count, datatype);
638 $mpi_scan(sendbuf, recvbuf, count, datatype, op, comm);
639 return 0;
640}
641
642int MPI_Exscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype,
643 MPI_Op op, MPI_Comm comm) {
644#ifndef _MPI_CONTRACT
645 $assert(_mpi_state == _MPI_INIT, "MPI_Exscan() cannot be invoked without "
646 "MPI_Init() being called before.\n");
647#endif
648 int place = $mpi_comm_place(comm);
649 int datatype_enum2int = (int)datatype;
650
651 // check consistency of a group of collective routine calls
652 $bundle checkerEntry = $mpi_create_coroutine_entry(EXSCAN_TAG, -1, -1, 1,
653 &datatype_enum2int);
654 $bundle specEntry = $mpi_check_collective(comm, checkerEntry);
655
656 $mpi_diff_coroutine_entries(specEntry, checkerEntry, place);
657 // MPI_IN_PLACE in this routine simply means that get sendbuf from recvbuf:
658 if (sendbuf != MPI_IN_PLACE)
659 $mpi_check_buffer(sendbuf, count, datatype);
660 if (place > 0)
661 $mpi_check_buffer(recvbuf, count, datatype);
662 $mpi_exscan(sendbuf, recvbuf, count, datatype, op, comm);
663 return 0;
664}
665
666/* ****************** End of collecitve routines ********************* */
667
668int MPI_Comm_dup(MPI_Comm comm, MPI_Comm * newcomm) {
669 $scope CMPI_PROC_SCOPE = $mpi_proc_scope(comm);
670
671#ifndef _MPI_CONTRACT
672 $assert(_mpi_state == _MPI_INIT, "MPI_Comm_dup() cannot be invoked without "
673 "MPI_Init() being called before.\n");
674#endif
675 $mpi_comm_dup(CMPI_PROC_SCOPE, comm, newcomm, "MPI_Comm_dup");
676 return 0;
677}
678
679int MPI_Comm_free(MPI_Comm * comm) {
680#ifndef _MPI_CONTRACT
681 $assert(_mpi_state == _MPI_INIT, "MPI_Comm_free() cannot be invoked without "
682 "MPI_Init() being called before.\n");
683#endif
684 $assert($is_derefable_pointer(comm), "The argument of MPI_Comm_free is NULL.");
685 $mpi_comm_free(comm, _mpi_state);
686 return 0;
687}
688
689int MPI_Init_thread( int *argc, char ***argv, int required, int *provided ){
690 _mpi_state = _MPI_INIT; //TODO: why set initialized flag here ??
691 *provided = MPI_THREAD_MULTIPLE;
692 return MPI_SUCCESS;
693}
694
695#endif
Note: See TracBrowser for help on using the repository browser.