Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Update copyright lines.
[simgrid.git] / src / smpi / bindings / smpi_f77_request.cpp
1 /* Copyright (c) 2010-2021. The SimGrid Team. All rights reserved.          */
2
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. */
5
6 #include "private.hpp"
7 #include "smpi_comm.hpp"
8 #include "smpi_datatype.hpp"
9 #include "smpi_request.hpp"
10
11
12 extern "C" { // This should really use the C linkage to be usable from Fortran
13
14 void mpi_send_init_(void *buf, int* count, int* datatype, int* dst, int* tag, int* comm, int* request, int* ierr) {
15   MPI_Request req;
16   *request = MPI_FORTRAN_REQUEST_NULL;
17   buf = static_cast<char *>(FORT_BOTTOM(buf));
18   *ierr = MPI_Send_init(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dst, *tag, simgrid::smpi::Comm::f2c(*comm), &req);
19   if(*ierr == MPI_SUCCESS && req != nullptr) {
20     *request = req->add_f();
21   }
22 }
23
24 void mpi_isend_(void *buf, int* count, int* datatype, int* dst, int* tag, int* comm, int* request, int* ierr) {
25   MPI_Request req;
26   *request = MPI_FORTRAN_REQUEST_NULL;
27   buf = static_cast<char *>(FORT_BOTTOM(buf));
28   *ierr = MPI_Isend(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dst, *tag, simgrid::smpi::Comm::f2c(*comm), &req);
29   if(*ierr == MPI_SUCCESS && req != nullptr) {
30     *request = req->add_f();
31   }
32 }
33
34 void mpi_irsend_(void *buf, int* count, int* datatype, int* dst, int* tag, int* comm, int* request, int* ierr) {
35   MPI_Request req;
36   *request = MPI_FORTRAN_REQUEST_NULL;
37   buf = static_cast<char *>(FORT_BOTTOM(buf));
38   *ierr = MPI_Irsend(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dst, *tag, simgrid::smpi::Comm::f2c(*comm), &req);
39   if(*ierr == MPI_SUCCESS && req != nullptr) {
40     *request = req->add_f();
41   }
42 }
43
44 void mpi_send_(void* buf, int* count, int* datatype, int* dst, int* tag, int* comm, int* ierr) {
45   buf = static_cast<char *>(FORT_BOTTOM(buf));
46    *ierr = MPI_Send(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dst, *tag, simgrid::smpi::Comm::f2c(*comm));
47 }
48
49 void mpi_rsend_(void* buf, int* count, int* datatype, int* dst, int* tag, int* comm, int* ierr) {
50   buf = static_cast<char *>(FORT_BOTTOM(buf));
51    *ierr = MPI_Rsend(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dst, *tag, simgrid::smpi::Comm::f2c(*comm));
52 }
53
54 void mpi_sendrecv_(void* sendbuf, int* sendcount, int* sendtype, int* dst, int* sendtag, void *recvbuf, int* recvcount,
55                    int* recvtype, int* src, int* recvtag, int* comm, MPI_Status* status, int* ierr) {
56   sendbuf = static_cast<char *>( FORT_BOTTOM(sendbuf));
57   recvbuf = static_cast<char *>( FORT_BOTTOM(recvbuf));
58    *ierr = MPI_Sendrecv(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype), *dst, *sendtag, recvbuf, *recvcount,
59                         simgrid::smpi::Datatype::f2c(*recvtype), *src, *recvtag, simgrid::smpi::Comm::f2c(*comm), FORT_STATUS_IGNORE(status));
60 }
61
62 void mpi_recv_init_(void *buf, int* count, int* datatype, int* src, int* tag, int* comm, int* request, int* ierr) {
63   MPI_Request req;
64   *request = MPI_FORTRAN_REQUEST_NULL;
65   buf = static_cast<char *>( FORT_BOTTOM(buf));
66   *ierr = MPI_Recv_init(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *src, *tag, simgrid::smpi::Comm::f2c(*comm), &req);
67   if(*ierr == MPI_SUCCESS) {
68     *request = req->add_f();
69   }
70 }
71
72 void mpi_irecv_(void *buf, int* count, int* datatype, int* src, int* tag, int* comm, int* request, int* ierr) {
73   MPI_Request req;
74   *request = MPI_FORTRAN_REQUEST_NULL;
75   buf = static_cast<char *>( FORT_BOTTOM(buf));
76   *ierr = MPI_Irecv(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *src, *tag, simgrid::smpi::Comm::f2c(*comm), &req);
77   if(*ierr == MPI_SUCCESS && req != nullptr) {
78     *request = req->add_f();
79   }
80 }
81
82 void mpi_recv_(void* buf, int* count, int* datatype, int* src, int* tag, int* comm, MPI_Status* status, int* ierr) {
83   buf = static_cast<char *>( FORT_BOTTOM(buf));
84   *ierr = MPI_Recv(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *src, *tag, simgrid::smpi::Comm::f2c(*comm), status);
85 }
86
87 void mpi_sendrecv_replace_ (void *buf, int* count, int* datatype, int* dst, int* sendtag, int* src, int* recvtag,
88                             int* comm, MPI_Status* status, int* ierr)
89 {
90   *ierr = MPI_Sendrecv_replace(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dst, *sendtag, *src,
91   *recvtag, simgrid::smpi::Comm::f2c(*comm), FORT_STATUS_IGNORE(status));
92 }
93
94 void mpi_ssend_ (void* buf, int* count, int* datatype, int* dest, int* tag, int* comm, int* ierr) {
95   *ierr = MPI_Ssend(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dest, *tag, simgrid::smpi::Comm::f2c(*comm));
96 }
97
98 void mpi_ssend_init_ (void* buf, int* count, int* datatype, int* dest, int* tag, int* comm, int* request, int* ierr) {
99   MPI_Request tmp;
100   *request = MPI_FORTRAN_REQUEST_NULL;
101   *ierr = MPI_Ssend_init(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dest, *tag, simgrid::smpi::Comm::f2c(*comm), &tmp);
102   if(*ierr == MPI_SUCCESS && tmp != nullptr) {
103     *request = tmp->add_f();
104   }
105 }
106
107 void mpi_bsend_ (void* buf, int* count, int* datatype, int *dest, int* tag, int* comm, int* ierr) {
108  *ierr = MPI_Bsend(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dest, *tag, simgrid::smpi::Comm::f2c(*comm));
109 }
110
111 void mpi_bsend_init_ (void* buf, int* count, int* datatype, int *dest, int* tag, int* comm, int*  request, int* ierr) {
112   MPI_Request tmp;
113   *request = MPI_FORTRAN_REQUEST_NULL;
114   *ierr = MPI_Bsend_init(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dest, *tag, simgrid::smpi::Comm::f2c(*comm), &tmp);
115   if(*ierr == MPI_SUCCESS && tmp != nullptr) {
116     *request = tmp->add_f();
117   }
118 }
119
120 void mpi_ibsend_ (void* buf, int* count, int* datatype, int *dest, int* tag, int* comm, int*  request, int* ierr) {
121   MPI_Request tmp;
122   *request = MPI_FORTRAN_REQUEST_NULL;
123   *ierr = MPI_Ibsend(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dest, *tag, simgrid::smpi::Comm::f2c(*comm), &tmp);
124   if(*ierr == MPI_SUCCESS && tmp != nullptr) {
125     *request = tmp->add_f();
126   }
127 }
128
129 void mpi_issend_ (void* buf, int* count, int* datatype, int *dest, int* tag, int* comm, int*  request, int* ierr) {
130   MPI_Request tmp;
131   *request = MPI_FORTRAN_REQUEST_NULL;
132   *ierr = MPI_Issend(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dest, *tag, simgrid::smpi::Comm::f2c(*comm), &tmp);
133   if(*ierr == MPI_SUCCESS && tmp != nullptr) {
134     *request = tmp->add_f();
135   }
136 }
137
138 void mpi_rsend_init_ (void* buf, int* count, int* datatype, int *dest, int* tag, int* comm, int*  request, int* ierr) {
139   MPI_Request tmp;
140   *request = MPI_FORTRAN_REQUEST_NULL;
141   *ierr = MPI_Rsend_init(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *dest, *tag, simgrid::smpi::Comm::f2c(*comm), &tmp);
142   if(*ierr == MPI_SUCCESS && tmp != nullptr) {
143     *request = tmp->add_f();
144   }
145 }
146
147 void mpi_start_(int* request, int* ierr) {
148   MPI_Request req = simgrid::smpi::Request::f2c(*request);
149
150   *ierr = MPI_Start(&req);
151 }
152
153 void mpi_startall_(int* count, int* requests, int* ierr)
154 {
155   std::vector<MPI_Request> reqs(*count);
156   for (int i = 0; i < *count; i++) {
157     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
158   }
159   *ierr = MPI_Startall(*count, reqs.data());
160 }
161
162 void mpi_wait_(int* request, MPI_Status* status, int* ierr) {
163    MPI_Request req = simgrid::smpi::Request::f2c(*request);
164
165    *ierr = MPI_Wait(&req, FORT_STATUS_IGNORE(status));
166    if(req==MPI_REQUEST_NULL){
167      simgrid::smpi::Request::free_f(*request);
168      *request=MPI_FORTRAN_REQUEST_NULL;
169    }
170 }
171
172 void mpi_waitany_(int* count, int* requests, int* index, MPI_Status* status, int* ierr)
173 {
174   std::vector<MPI_Request> reqs(*count);
175   for (int i = 0; i < *count; i++) {
176     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
177   }
178   *ierr = MPI_Waitany(*count, reqs.data(), index, status);
179   if(*index!=MPI_UNDEFINED){
180     if(reqs[*index]==MPI_REQUEST_NULL){
181         simgrid::smpi::Request::free_f(requests[*index]);
182         requests[*index]=MPI_FORTRAN_REQUEST_NULL;
183     }
184   *index=*index+1;
185   }
186 }
187
188 void mpi_waitall_(int* count, int* requests, MPI_Status* status, int* ierr)
189 {
190   std::vector<MPI_Request> reqs(*count);
191   for (int i = 0; i < *count; i++) {
192     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
193   }
194   *ierr = MPI_Waitall(*count, reqs.data(), FORT_STATUSES_IGNORE(status));
195   for (int i = 0; i < *count; i++) {
196     if (reqs[i] == MPI_REQUEST_NULL) {
197       simgrid::smpi::Request::free_f(requests[i]);
198       requests[i] = MPI_FORTRAN_REQUEST_NULL;
199     }
200   }
201 }
202
203 void mpi_waitsome_ (int* incount, int* requests, int *outcount, int *indices, MPI_Status* status, int* ierr)
204 {
205   std::vector<MPI_Request> reqs(*incount);
206   for (int i = 0; i < *incount; i++) {
207     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
208   }
209   *ierr = MPI_Waitsome(*incount, reqs.data(), outcount, indices, status);
210   for (int 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;
214     }
215     indices[i]++;
216   }
217 }
218
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;
225   }
226 }
227
228 void mpi_testall_(int* count, int* requests, int* flag, MPI_Status* statuses, int* ierr)
229 {
230   std::vector<MPI_Request> reqs(*count);
231   for (int i = 0; i < *count; i++) {
232     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
233   }
234   *ierr = MPI_Testall(*count, reqs.data(), flag, FORT_STATUSES_IGNORE(statuses));
235   for (int 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;
239     }
240   }
241 }
242
243 void mpi_testany_ (int* count, int* requests, int *index, int *flag, MPI_Status* status, int* ierr)
244 {
245   std::vector<MPI_Request> reqs(*count);
246   for (int i = 0; i < *count; i++) {
247     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
248   }
249   *ierr = MPI_Testany(*count, reqs.data(), index, flag, FORT_STATUS_IGNORE(status));
250   if(*index!=MPI_UNDEFINED){
251     if(reqs[*index]==MPI_REQUEST_NULL){
252     simgrid::smpi::Request::free_f(requests[*index]);
253     requests[*index]=MPI_FORTRAN_REQUEST_NULL;
254     }
255   *index=*index+1;
256   }
257 }
258
259 void mpi_testsome_(int* incount, int* requests, int* outcount, int* indices, MPI_Status* statuses, int* ierr)
260 {
261   std::vector<MPI_Request> reqs(*incount);
262   for (int i = 0; i < *incount; i++) {
263     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
264     indices[i]=0;
265   }
266   *ierr = MPI_Testsome(*incount, reqs.data(), outcount, indices, FORT_STATUSES_IGNORE(statuses));
267   for (int i = 0; i < *incount; i++) {
268     if(reqs[indices[i]]==MPI_REQUEST_NULL){
269       simgrid::smpi::Request::free_f(requests[indices[i]]);
270       requests[indices[i]]=MPI_FORTRAN_REQUEST_NULL;
271     }
272     indices[i]++;
273   }
274 }
275
276 void mpi_probe_ (int* source, int* tag, int* comm, MPI_Status*  status, int* ierr) {
277  *ierr = MPI_Probe(*source, *tag, simgrid::smpi::Comm::f2c(*comm), FORT_STATUS_IGNORE(status));
278 }
279
280
281 void mpi_iprobe_ (int* source, int* tag, int* comm, int* flag, MPI_Status*  status, int* ierr) {
282  *ierr = MPI_Iprobe(*source, *tag, simgrid::smpi::Comm::f2c(*comm), flag, status);
283 }
284
285 }