Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
More uses of std::make_unique.
[simgrid.git] / src / smpi / bindings / smpi_f77_request.cpp
1 /* Copyright (c) 2010-2020. 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   MPI_Request* reqs;
155
156   reqs = xbt_new(MPI_Request, *count);
157   for (int i = 0; i < *count; i++) {
158     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
159   }
160   *ierr = MPI_Startall(*count, reqs);
161   xbt_free(reqs);
162 }
163
164 void mpi_wait_(int* request, MPI_Status* status, int* ierr) {
165    MPI_Request req = simgrid::smpi::Request::f2c(*request);
166
167    *ierr = MPI_Wait(&req, FORT_STATUS_IGNORE(status));
168    if(req==MPI_REQUEST_NULL){
169      simgrid::smpi::Request::free_f(*request);
170      *request=MPI_FORTRAN_REQUEST_NULL;
171    }
172 }
173
174 void mpi_waitany_(int* count, int* requests, int* index, MPI_Status* status, int* ierr) {
175   MPI_Request* reqs;
176
177   reqs = xbt_new(MPI_Request, *count);
178   for (int i = 0; i < *count; i++) {
179     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
180   }
181   *ierr = MPI_Waitany(*count, reqs, index, status);
182   if(*index!=MPI_UNDEFINED){
183     if(reqs[*index]==MPI_REQUEST_NULL){
184         simgrid::smpi::Request::free_f(requests[*index]);
185         requests[*index]=MPI_FORTRAN_REQUEST_NULL;
186     }
187   *index=*index+1;
188   }
189   xbt_free(reqs);
190 }
191
192 void mpi_waitall_(int* count, int* requests, MPI_Status* status, int* ierr) {
193   MPI_Request* reqs;
194   int i;
195
196   reqs = xbt_new(MPI_Request, *count);
197   for(i = 0; i < *count; i++) {
198     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
199   }
200   *ierr = MPI_Waitall(*count, reqs, FORT_STATUSES_IGNORE(status));
201   for(i = 0; i < *count; i++) {
202       if(reqs[i]==MPI_REQUEST_NULL){
203           simgrid::smpi::Request::free_f(requests[i]);
204           requests[i]=MPI_FORTRAN_REQUEST_NULL;
205       }
206   }
207
208   xbt_free(reqs);
209 }
210
211 void mpi_waitsome_ (int* incount, int* requests, int *outcount, int *indices, MPI_Status* status, int* ierr)
212 {
213   MPI_Request* reqs;
214   int i;
215
216   reqs = xbt_new(MPI_Request, *incount);
217   for(i = 0; i < *incount; i++) {
218     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
219   }
220   *ierr = MPI_Waitsome(*incount, reqs, outcount, indices, status);
221   for(i=0;i<*outcount;i++){
222     if(reqs[indices[i]]==MPI_REQUEST_NULL){
223         simgrid::smpi::Request::free_f(requests[indices[i]]);
224         requests[indices[i]]=MPI_FORTRAN_REQUEST_NULL;
225     }
226     indices[i]++;
227   }
228   xbt_free(reqs);
229 }
230
231 void mpi_test_ (int * request, int *flag, MPI_Status * status, int* ierr){
232   MPI_Request req = simgrid::smpi::Request::f2c(*request);
233   *ierr= MPI_Test(&req, flag, FORT_STATUS_IGNORE(status));
234   if(req==MPI_REQUEST_NULL){
235       simgrid::smpi::Request::free_f(*request);
236       *request=MPI_FORTRAN_REQUEST_NULL;
237   }
238 }
239
240 void mpi_testall_ (int* count, int * requests,  int *flag, MPI_Status * statuses, int* ierr){
241   int i;
242   MPI_Request* reqs = xbt_new(MPI_Request, *count);
243   for(i = 0; i < *count; i++) {
244     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
245   }
246   *ierr= MPI_Testall(*count, reqs, flag, FORT_STATUSES_IGNORE(statuses));
247   for(i = 0; i < *count; i++) {
248     if(reqs[i]==MPI_REQUEST_NULL){
249         simgrid::smpi::Request::free_f(requests[i]);
250         requests[i]=MPI_FORTRAN_REQUEST_NULL;
251     }
252   }
253   xbt_free(reqs);
254 }
255
256 void mpi_testany_ (int* count, int* requests, int *index, int *flag, MPI_Status* status, int* ierr)
257 {
258   MPI_Request* reqs;
259
260   reqs = xbt_new(MPI_Request, *count);
261   for (int i = 0; i < *count; i++) {
262     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
263   }
264   *ierr = MPI_Testany(*count, reqs, index, flag, FORT_STATUS_IGNORE(status));
265   if(*index!=MPI_UNDEFINED){
266     if(reqs[*index]==MPI_REQUEST_NULL){
267     simgrid::smpi::Request::free_f(requests[*index]);
268     requests[*index]=MPI_FORTRAN_REQUEST_NULL;
269     }
270   *index=*index+1;
271   }
272   xbt_free(reqs);
273 }
274
275 void mpi_testsome_ (int* incount, int*  requests, int* outcount, int* indices, MPI_Status*  statuses, int* ierr) {
276   MPI_Request* reqs;
277   int i;
278
279   reqs = xbt_new(MPI_Request, *incount);
280   for(i = 0; i < *incount; i++) {
281     reqs[i] = simgrid::smpi::Request::f2c(requests[i]);
282     indices[i]=0;
283   }
284   *ierr = MPI_Testsome(*incount, reqs, outcount, indices, FORT_STATUSES_IGNORE(statuses));
285   for(i=0;i<*incount;i++){
286     if(reqs[indices[i]]==MPI_REQUEST_NULL){
287       simgrid::smpi::Request::free_f(requests[indices[i]]);
288       requests[indices[i]]=MPI_FORTRAN_REQUEST_NULL;
289     }
290     indices[i]++;
291   }
292   xbt_free(reqs);
293 }
294
295 void mpi_probe_ (int* source, int* tag, int* comm, MPI_Status*  status, int* ierr) {
296  *ierr = MPI_Probe(*source, *tag, simgrid::smpi::Comm::f2c(*comm), FORT_STATUS_IGNORE(status));
297 }
298
299
300 void mpi_iprobe_ (int* source, int* tag, int* comm, int* flag, MPI_Status*  status, int* ierr) {
301  *ierr = MPI_Iprobe(*source, *tag, simgrid::smpi::Comm::f2c(*comm), flag, status);
302 }
303
304 }