1 /* Copyright (c) 2010-2018. The SimGrid Team. All rights reserved. */
3 /* This program is free software; you can redistribute it and/or modify it
4 * under the terms of the license (GNU LGPL) which comes with this package. */
7 #include "smpi_comm.hpp"
8 #include "smpi_datatype.hpp"
9 #include "smpi_request.hpp"
12 extern "C" { // This should really use the C linkage to be usable from Fortran
14 void mpi_send_init_(void *buf, int* count, int* datatype, int* dst, int* tag, int* comm, int* request, int* ierr) {
16 buf = static_cast<char *>(FORT_BOTTOM(buf));
17 *ierr = MPI_Send_init(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dst, *tag, simgrid::smpi::Comm::f2c(*comm), &req);
18 if(*ierr == MPI_SUCCESS) {
19 *request = req->add_f();
23 void mpi_isend_(void *buf, int* count, int* datatype, int* dst, int* tag, int* comm, int* request, int* ierr) {
25 buf = static_cast<char *>(FORT_BOTTOM(buf));
26 *ierr = MPI_Isend(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dst, *tag, simgrid::smpi::Comm::f2c(*comm), &req);
27 if(*ierr == MPI_SUCCESS) {
28 *request = req->add_f();
32 void mpi_irsend_(void *buf, int* count, int* datatype, int* dst, int* tag, int* comm, int* request, int* ierr) {
34 buf = static_cast<char *>(FORT_BOTTOM(buf));
35 *ierr = MPI_Irsend(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dst, *tag, simgrid::smpi::Comm::f2c(*comm), &req);
36 if(*ierr == MPI_SUCCESS) {
37 *request = req->add_f();
41 void mpi_send_(void* buf, int* count, int* datatype, int* dst, int* tag, int* comm, int* ierr) {
42 buf = static_cast<char *>(FORT_BOTTOM(buf));
43 *ierr = MPI_Send(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dst, *tag, simgrid::smpi::Comm::f2c(*comm));
46 void mpi_rsend_(void* buf, int* count, int* datatype, int* dst, int* tag, int* comm, int* ierr) {
47 buf = static_cast<char *>(FORT_BOTTOM(buf));
48 *ierr = MPI_Rsend(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dst, *tag, simgrid::smpi::Comm::f2c(*comm));
51 void mpi_sendrecv_(void* sendbuf, int* sendcount, int* sendtype, int* dst, int* sendtag, void *recvbuf, int* recvcount,
52 int* recvtype, int* src, int* recvtag, int* comm, MPI_Status* status, int* ierr) {
53 sendbuf = static_cast<char *>( FORT_BOTTOM(sendbuf));
54 recvbuf = static_cast<char *>( FORT_BOTTOM(recvbuf));
55 *ierr = MPI_Sendrecv(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype), *dst, *sendtag, recvbuf, *recvcount,
56 simgrid::smpi::Datatype::f2c(*recvtype), *src, *recvtag, simgrid::smpi::Comm::f2c(*comm), FORT_STATUS_IGNORE(status));
59 void mpi_recv_init_(void *buf, int* count, int* datatype, int* src, int* tag, int* comm, int* request, int* ierr) {
61 buf = static_cast<char *>( FORT_BOTTOM(buf));
62 *ierr = MPI_Recv_init(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *src, *tag, simgrid::smpi::Comm::f2c(*comm), &req);
63 if(*ierr == MPI_SUCCESS) {
64 *request = req->add_f();
68 void mpi_irecv_(void *buf, int* count, int* datatype, int* src, int* tag, int* comm, int* request, int* ierr) {
70 buf = static_cast<char *>( FORT_BOTTOM(buf));
71 *ierr = MPI_Irecv(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *src, *tag, simgrid::smpi::Comm::f2c(*comm), &req);
72 if(*ierr == MPI_SUCCESS) {
73 *request = req->add_f();
77 void mpi_recv_(void* buf, int* count, int* datatype, int* src, int* tag, int* comm, MPI_Status* status, int* ierr) {
78 buf = static_cast<char *>( FORT_BOTTOM(buf));
79 *ierr = MPI_Recv(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *src, *tag, simgrid::smpi::Comm::f2c(*comm), status);
82 void mpi_sendrecv_replace_ (void *buf, int* count, int* datatype, int* dst, int* sendtag, int* src, int* recvtag,
83 int* comm, MPI_Status* status, int* ierr)
85 *ierr = MPI_Sendrecv_replace(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dst, *sendtag, *src,
86 *recvtag, simgrid::smpi::Comm::f2c(*comm), FORT_STATUS_IGNORE(status));
89 void mpi_ssend_ (void* buf, int* count, int* datatype, int* dest, int* tag, int* comm, int* ierr) {
90 *ierr = MPI_Ssend(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dest, *tag, simgrid::smpi::Comm::f2c(*comm));
93 void mpi_ssend_init_ (void* buf, int* count, int* datatype, int* dest, int* tag, int* comm, int* request, int* ierr) {
95 *ierr = MPI_Ssend_init(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dest, *tag, simgrid::smpi::Comm::f2c(*comm), &tmp);
96 if(*ierr == MPI_SUCCESS) {
97 *request = tmp->add_f();
101 void mpi_bsend_ (void* buf, int* count, int* datatype, int *dest, int* tag, int* comm, int* ierr) {
102 *ierr = MPI_Bsend(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dest, *tag, simgrid::smpi::Comm::f2c(*comm));
105 void mpi_bsend_init_ (void* buf, int* count, int* datatype, int *dest, int* tag, int* comm, int* request, int* ierr) {
107 *ierr = MPI_Bsend_init(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dest, *tag, simgrid::smpi::Comm::f2c(*comm), &tmp);
108 if(*ierr == MPI_SUCCESS) {
109 *request = tmp->add_f();
113 void mpi_ibsend_ (void* buf, int* count, int* datatype, int *dest, int* tag, int* comm, int* request, int* ierr) {
115 *ierr = MPI_Ibsend(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dest, *tag, simgrid::smpi::Comm::f2c(*comm), &tmp);
116 if(*ierr == MPI_SUCCESS) {
117 *request = tmp->add_f();
121 void mpi_issend_ (void* buf, int* count, int* datatype, int *dest, int* tag, int* comm, int* request, int* ierr) {
123 *ierr = MPI_Issend(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dest, *tag, simgrid::smpi::Comm::f2c(*comm), &tmp);
124 if(*ierr == MPI_SUCCESS) {
125 *request = tmp->add_f();
129 void mpi_rsend_init_ (void* buf, int* count, int* datatype, int *dest, int* tag, int* comm, int* request, int* ierr) {
131 *ierr = MPI_Rsend_init(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dest, *tag, simgrid::smpi::Comm::f2c(*comm), &tmp);
132 if(*ierr == MPI_SUCCESS) {
133 *request = tmp->add_f();
137 void mpi_start_(int* request, int* ierr) {
138 MPI_Request req = simgrid::smpi::Request::f2c(*request);
140 *ierr = MPI_Start(&req);
143 void mpi_startall_(int* count, int* requests, int* ierr) {
147 reqs = xbt_new(MPI_Request, *count);
148 for(i = 0; i < *count; i++) {
149 reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
151 *ierr = MPI_Startall(*count, reqs);
155 void mpi_wait_(int* request, MPI_Status* status, int* ierr) {
156 MPI_Request req = simgrid::smpi::Request::f2c(*request);
158 *ierr = MPI_Wait(&req, FORT_STATUS_IGNORE(status));
159 if(req==MPI_REQUEST_NULL){
160 simgrid::smpi::Request::free_f(*request);
161 *request=MPI_FORTRAN_REQUEST_NULL;
165 void mpi_waitany_(int* count, int* requests, int* index, MPI_Status* status, int* ierr) {
169 reqs = xbt_new(MPI_Request, *count);
170 for(i = 0; i < *count; i++) {
171 reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
173 *ierr = MPI_Waitany(*count, reqs, index, status);
174 if(reqs[*index]==MPI_REQUEST_NULL){
175 simgrid::smpi::Request::free_f(requests[*index]);
176 requests[*index]=MPI_FORTRAN_REQUEST_NULL;
181 void mpi_waitall_(int* count, int* requests, MPI_Status* status, int* ierr) {
185 reqs = xbt_new(MPI_Request, *count);
186 for(i = 0; i < *count; i++) {
187 reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
189 *ierr = MPI_Waitall(*count, reqs, FORT_STATUSES_IGNORE(status));
190 for(i = 0; i < *count; i++) {
191 if(reqs[i]==MPI_REQUEST_NULL){
192 simgrid::smpi::Request::free_f(requests[i]);
193 requests[i]=MPI_FORTRAN_REQUEST_NULL;
200 void mpi_waitsome_ (int* incount, int* requests, int *outcount, int *indices, MPI_Status* status, int* ierr)
205 reqs = xbt_new(MPI_Request, *incount);
206 for(i = 0; i < *incount; i++) {
207 reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
209 *ierr = MPI_Waitsome(*incount, reqs, outcount, indices, status);
210 for(i=0;i<*outcount;i++){
211 if(reqs[indices[i]]==MPI_REQUEST_NULL){
212 simgrid::smpi::Request::free_f(requests[indices[i]]);
213 requests[indices[i]]=MPI_FORTRAN_REQUEST_NULL;
219 void mpi_test_ (int * request, int *flag, MPI_Status * status, int* ierr){
220 MPI_Request req = simgrid::smpi::Request::f2c(*request);
221 *ierr= MPI_Test(&req, flag, FORT_STATUS_IGNORE(status));
222 if(req==MPI_REQUEST_NULL){
223 simgrid::smpi::Request::free_f(*request);
224 *request=MPI_FORTRAN_REQUEST_NULL;
228 void mpi_testall_ (int* count, int * requests, int *flag, MPI_Status * statuses, int* ierr){
230 MPI_Request* reqs = xbt_new(MPI_Request, *count);
231 for(i = 0; i < *count; i++) {
232 reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
234 *ierr= MPI_Testall(*count, reqs, flag, FORT_STATUSES_IGNORE(statuses));
235 for(i = 0; i < *count; i++) {
236 if(reqs[i]==MPI_REQUEST_NULL){
237 simgrid::smpi::Request::free_f(requests[i]);
238 requests[i]=MPI_FORTRAN_REQUEST_NULL;
244 void mpi_testany_ (int* count, int* requests, int *index, int *flag, MPI_Status* status, int* ierr)
249 reqs = xbt_new(MPI_Request, *count);
250 for(i = 0; i < *count; i++) {
251 reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
253 *ierr = MPI_Testany(*count, reqs, index, flag, FORT_STATUS_IGNORE(status));
254 if(*index!=MPI_UNDEFINED && reqs[*index]==MPI_REQUEST_NULL){
255 simgrid::smpi::Request::free_f(requests[*index]);
256 requests[*index]=MPI_FORTRAN_REQUEST_NULL;
261 void mpi_testsome_ (int* incount, int* requests, int* outcount, int* indices, MPI_Status* statuses, int* ierr) {
265 reqs = xbt_new(MPI_Request, *incount);
266 for(i = 0; i < *incount; i++) {
267 reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
270 *ierr = MPI_Testsome(*incount, reqs, outcount, indices, FORT_STATUSES_IGNORE(statuses));
271 for(i=0;i<*incount;i++){
272 if(indices[i] && reqs[indices[i]]==MPI_REQUEST_NULL){
273 simgrid::smpi::Request::free_f(requests[indices[i]]);
274 requests[indices[i]]=MPI_FORTRAN_REQUEST_NULL;
280 void mpi_probe_ (int* source, int* tag, int* comm, MPI_Status* status, int* ierr) {
281 *ierr = MPI_Probe(*source, *tag, simgrid::smpi::Comm::f2c(*comm), FORT_STATUS_IGNORE(status));
285 void mpi_iprobe_ (int* source, int* tag, int* comm, int* flag, MPI_Status* status, int* ierr) {
286 *ierr = MPI_Iprobe(*source, *tag, simgrid::smpi::Comm::f2c(*comm), flag, status);