Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
50296432f614a7daf114dcdf0c528fee0d225c28
[simgrid.git] / src / smpi / bindings / smpi_f77_coll.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_coll.hpp"
8 #include "smpi_comm.hpp"
9 #include "smpi_datatype.hpp"
10 #include "smpi_op.hpp"
11 #include "smpi_request.hpp"
12
13 extern "C" { // This should really use the C linkage to be usable from Fortran
14
15 void mpi_barrier_(int* comm, int* ierr) {
16   *ierr = MPI_Barrier(simgrid::smpi::Comm::f2c(*comm));
17 }
18
19 void mpi_bcast_(void *buf, int* count, int* datatype, int* root, int* comm, int* ierr) {
20   *ierr = MPI_Bcast(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *root, simgrid::smpi::Comm::f2c(*comm));
21 }
22
23 void mpi_reduce_(void* sendbuf, void* recvbuf, int* count, int* datatype, int* op, int* root, int* comm, int* ierr) {
24   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
25   sendbuf = static_cast<char *>( FORT_BOTTOM(sendbuf));
26   recvbuf = static_cast<char *>( FORT_BOTTOM(recvbuf));
27   *ierr = MPI_Reduce(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), *root, simgrid::smpi::Comm::f2c(*comm));
28 }
29
30 void mpi_allreduce_(void* sendbuf, void* recvbuf, int* count, int* datatype, int* op, int* comm, int* ierr) {
31   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
32   *ierr = MPI_Allreduce(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm));
33 }
34
35 void mpi_reduce_scatter_(void* sendbuf, void* recvbuf, int* recvcounts, int* datatype, int* op, int* comm, int* ierr) {
36   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
37   *ierr = MPI_Reduce_scatter(sendbuf, recvbuf, recvcounts, simgrid::smpi::Datatype::f2c(*datatype),
38                         simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm));
39 }
40
41 void mpi_scatter_(void* sendbuf, int* sendcount, int* sendtype, void* recvbuf, int* recvcount, int* recvtype,
42                    int* root, int* comm, int* ierr) {
43   recvbuf = static_cast<char *>( FORT_IN_PLACE(recvbuf));
44   *ierr = MPI_Scatter(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
45                       recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm));
46 }
47
48 void mpi_scatterv_(void* sendbuf, int* sendcounts, int* displs, int* sendtype,
49                    void* recvbuf, int* recvcount, int* recvtype, int* root, int* comm, int* ierr) {
50   recvbuf = static_cast<char *>( FORT_IN_PLACE(recvbuf));
51   *ierr = MPI_Scatterv(sendbuf, sendcounts, displs, simgrid::smpi::Datatype::f2c(*sendtype),
52                       recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm));
53 }
54
55 void mpi_gather_(void* sendbuf, int* sendcount, int* sendtype, void* recvbuf, int* recvcount, int* recvtype,
56                   int* root, int* comm, int* ierr) {
57   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
58   sendbuf = sendbuf!=MPI_IN_PLACE ? static_cast<char *>( FORT_BOTTOM(sendbuf)) : MPI_IN_PLACE;
59   recvbuf = static_cast<char *>( FORT_BOTTOM(recvbuf));
60   *ierr = MPI_Gather(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
61                      recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm));
62 }
63
64 void mpi_gatherv_(void* sendbuf, int* sendcount, int* sendtype,
65                   void* recvbuf, int* recvcounts, int* displs, int* recvtype, int* root, int* comm, int* ierr) {
66   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
67   sendbuf = sendbuf!=MPI_IN_PLACE ? static_cast<char *>( FORT_BOTTOM(sendbuf)) : MPI_IN_PLACE;
68   recvbuf = static_cast<char *>( FORT_BOTTOM(recvbuf));
69   *ierr = MPI_Gatherv(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
70                      recvbuf, recvcounts, displs, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm));
71 }
72
73 void mpi_allgather_(void* sendbuf, int* sendcount, int* sendtype, void* recvbuf, int* recvcount, int* recvtype,
74                      int* comm, int* ierr) {
75   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
76   *ierr = MPI_Allgather(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
77                         recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm));
78 }
79
80 void mpi_allgatherv_(void* sendbuf, int* sendcount, int* sendtype,
81                      void* recvbuf, int* recvcounts,int* displs, int* recvtype, int* comm, int* ierr) {
82   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
83   *ierr = MPI_Allgatherv(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
84                         recvbuf, recvcounts, displs, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm));
85 }
86
87 void mpi_scan_(void* sendbuf, void* recvbuf, int* count, int* datatype, int* op, int* comm, int* ierr) {
88   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
89   *ierr = MPI_Scan(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype),
90                    simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm));
91 }
92
93 void mpi_alltoall_(void* sendbuf, int* sendcount, int* sendtype,
94                     void* recvbuf, int* recvcount, int* recvtype, int* comm, int* ierr) {
95   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
96   *ierr = MPI_Alltoall(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
97                        recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm));
98 }
99
100 void mpi_alltoallv_(void* sendbuf, int* sendcounts, int* senddisps, int* sendtype,
101                     void* recvbuf, int* recvcounts, int* recvdisps, int* recvtype, int* comm, int* ierr) {
102   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
103   *ierr = MPI_Alltoallv(sendbuf, sendcounts, senddisps, simgrid::smpi::Datatype::f2c(*sendtype),
104                        recvbuf, recvcounts, recvdisps, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm));
105 }
106
107 void mpi_reduce_local_ (void *inbuf, void *inoutbuf, int* count, int* datatype, int* op, int* ierr){
108  *ierr = MPI_Reduce_local(inbuf, inoutbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op));
109 }
110
111 void mpi_reduce_scatter_block_ (void *sendbuf, void *recvbuf, int* recvcount, int* datatype, int* op, int* comm,
112                                 int* ierr)
113 {
114   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
115  *ierr = MPI_Reduce_scatter_block(sendbuf, recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op),
116                                   simgrid::smpi::Comm::f2c(*comm));
117 }
118
119 void mpi_alltoallw_ ( void *sendbuf, int *sendcnts, int *sdispls, int* old_sendtypes, void *recvbuf, int *recvcnts,
120                       int *rdispls, int* old_recvtypes, int* comm, int* ierr){
121   int size = simgrid::smpi::Comm::f2c(*comm)->size();
122   MPI_Datatype* sendtypes = new MPI_Datatype[size];
123   MPI_Datatype* recvtypes = new MPI_Datatype[size];
124   for(int i=0; i< size; i++){
125     if(FORT_IN_PLACE(sendbuf)!=MPI_IN_PLACE)
126       sendtypes[i] = simgrid::smpi::Datatype::f2c(old_sendtypes[i]);
127     recvtypes[i] = simgrid::smpi::Datatype::f2c(old_recvtypes[i]);
128   }
129   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
130  *ierr = MPI_Alltoallw( sendbuf, sendcnts, sdispls, sendtypes, recvbuf, recvcnts, rdispls,
131                         recvtypes, simgrid::smpi::Comm::f2c(*comm));
132   delete[] sendtypes;
133   delete[] recvtypes;
134 }
135
136 void mpi_exscan_ (void *sendbuf, void *recvbuf, int* count, int* datatype, int* op, int* comm, int* ierr){
137  *ierr = MPI_Exscan(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm));
138 }
139
140 void mpi_ibarrier_(int* comm, int* request, int* ierr) {
141   MPI_Request req;
142   *ierr = MPI_Ibarrier(simgrid::smpi::Comm::f2c(*comm), &req);
143   if(*ierr == MPI_SUCCESS) {
144     *request = req->add_f();
145   }
146 }
147
148 void mpi_ibcast_(void *buf, int* count, int* datatype, int* root, int* comm, int* request, int* ierr) {
149   MPI_Request req;
150   *ierr = MPI_Ibcast(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *root, simgrid::smpi::Comm::f2c(*comm), &req);
151   if(*ierr == MPI_SUCCESS) {
152     *request = req->add_f();
153   }
154 }
155
156 void mpi_ireduce_(void* sendbuf, void* recvbuf, int* count, int* datatype, int* op, int* root, int* comm, int* request, int* ierr) {
157   MPI_Request req;
158   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
159   sendbuf = static_cast<char *>( FORT_BOTTOM(sendbuf));
160   recvbuf = static_cast<char *>( FORT_BOTTOM(recvbuf));
161   *ierr = MPI_Ireduce(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), *root, simgrid::smpi::Comm::f2c(*comm), &req);
162   if(*ierr == MPI_SUCCESS) {
163     *request = req->add_f();
164   }
165 }
166
167 void mpi_iallreduce_(void* sendbuf, void* recvbuf, int* count, int* datatype, int* op, int* comm, int* request, int* ierr) {
168   MPI_Request req;
169   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
170   *ierr = MPI_Iallreduce(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm), &req);
171   if(*ierr == MPI_SUCCESS) {
172     *request = req->add_f();
173   }
174 }
175
176 void mpi_ireduce_scatter_(void* sendbuf, void* recvbuf, int* recvcounts, int* datatype, int* op, int* comm, int* request, int* ierr) {
177   MPI_Request req;
178   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
179   *ierr = MPI_Ireduce_scatter(sendbuf, recvbuf, recvcounts, simgrid::smpi::Datatype::f2c(*datatype),
180                         simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm), &req);
181   if(*ierr == MPI_SUCCESS) {
182     *request = req->add_f();
183   }
184 }
185
186 void mpi_iscatter_(void* sendbuf, int* sendcount, int* sendtype, void* recvbuf, int* recvcount, int* recvtype,
187                    int* root, int* comm, int* request, int* ierr) {
188   MPI_Request req;
189   recvbuf = static_cast<char *>( FORT_IN_PLACE(recvbuf));
190   *ierr = MPI_Iscatter(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
191                       recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm), &req);
192   if(*ierr == MPI_SUCCESS) {
193     *request = req->add_f();
194   }
195 }
196
197 void mpi_iscatterv_(void* sendbuf, int* sendcounts, int* displs, int* sendtype,
198                    void* recvbuf, int* recvcount, int* recvtype, int* root, int* comm, int* request, int* ierr) {
199   MPI_Request req;
200   recvbuf = static_cast<char *>( FORT_IN_PLACE(recvbuf));
201   *ierr = MPI_Iscatterv(sendbuf, sendcounts, displs, simgrid::smpi::Datatype::f2c(*sendtype),
202                       recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm), &req);
203   if(*ierr == MPI_SUCCESS) {
204     *request = req->add_f();
205   }
206 }
207
208 void mpi_igather_(void* sendbuf, int* sendcount, int* sendtype, void* recvbuf, int* recvcount, int* recvtype,
209                   int* root, int* comm, int* request, int* ierr) {
210   MPI_Request req;
211   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
212   sendbuf = sendbuf!=MPI_IN_PLACE ? static_cast<char *>( FORT_BOTTOM(sendbuf)) : MPI_IN_PLACE;
213   recvbuf = static_cast<char *>( FORT_BOTTOM(recvbuf));
214   *ierr = MPI_Igather(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
215                      recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm), &req);
216   if(*ierr == MPI_SUCCESS) {
217     *request = req->add_f();
218   }
219 }
220
221 void mpi_igatherv_(void* sendbuf, int* sendcount, int* sendtype,
222                   void* recvbuf, int* recvcounts, int* displs, int* recvtype, int* root, int* comm, int* request, int* ierr) {
223   MPI_Request req;
224   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
225   sendbuf = sendbuf!=MPI_IN_PLACE ? static_cast<char *>( FORT_BOTTOM(sendbuf)) : MPI_IN_PLACE;
226   recvbuf = static_cast<char *>( FORT_BOTTOM(recvbuf));
227   *ierr = MPI_Igatherv(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
228                      recvbuf, recvcounts, displs, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm), &req);
229   if(*ierr == MPI_SUCCESS) {
230     *request = req->add_f();
231   }
232 }
233
234 void mpi_iallgather_(void* sendbuf, int* sendcount, int* sendtype, void* recvbuf, int* recvcount, int* recvtype,
235                      int* comm, int* request, int* ierr) {
236   MPI_Request req;
237   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
238   *ierr = MPI_Iallgather(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
239                         recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm), &req);
240   if(*ierr == MPI_SUCCESS) {
241     *request = req->add_f();
242   }
243 }
244
245 void mpi_iallgatherv_(void* sendbuf, int* sendcount, int* sendtype,
246                      void* recvbuf, int* recvcounts,int* displs, int* recvtype, int* comm, int* request, int* ierr) {
247   MPI_Request req;
248   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
249   *ierr = MPI_Iallgatherv(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
250                         recvbuf, recvcounts, displs, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm), &req);
251   if(*ierr == MPI_SUCCESS) {
252     *request = req->add_f();
253   }
254 }
255
256 void mpi_iscan_(void* sendbuf, void* recvbuf, int* count, int* datatype, int* op, int* comm, int* request, int* ierr) {
257   MPI_Request req;
258   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
259   *ierr = MPI_Iscan(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype),
260                    simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm), &req);
261   if(*ierr == MPI_SUCCESS) {
262     *request = req->add_f();
263   }
264 }
265
266 void mpi_ialltoall_(void* sendbuf, int* sendcount, int* sendtype,
267                     void* recvbuf, int* recvcount, int* recvtype, int* comm, int* request, int* ierr) {
268   MPI_Request req;
269   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
270   *ierr = MPI_Ialltoall(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
271                        recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm), &req);
272   if(*ierr == MPI_SUCCESS) {
273     *request = req->add_f();
274   }
275 }
276
277 void mpi_ialltoallv_(void* sendbuf, int* sendcounts, int* senddisps, int* sendtype,
278                     void* recvbuf, int* recvcounts, int* recvdisps, int* recvtype, int* comm, int* request, int* ierr) {
279   MPI_Request req;
280   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
281   *ierr = MPI_Ialltoallv(sendbuf, sendcounts, senddisps, simgrid::smpi::Datatype::f2c(*sendtype),
282                        recvbuf, recvcounts, recvdisps, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm), &req);
283   if(*ierr == MPI_SUCCESS) {
284     *request = req->add_f();
285   }
286 }
287
288 void mpi_ireduce_scatter_block_ (void *sendbuf, void *recvbuf, int* recvcount, int* datatype, int* op, int* comm,
289                                  int* request, int* ierr)
290 {
291   MPI_Request req;
292   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
293  *ierr = MPI_Ireduce_scatter_block(sendbuf, recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op),
294                                   simgrid::smpi::Comm::f2c(*comm), &req);
295   if(*ierr == MPI_SUCCESS) {
296     *request = req->add_f();
297   }
298 }
299
300 void mpi_ialltoallw_ ( void *sendbuf, int *sendcnts, int *sdispls, int* old_sendtypes, void *recvbuf, int *recvcnts,
301                       int *rdispls, int* old_recvtypes, int* comm, int* request, int* ierr){
302   MPI_Request req;
303   int size = simgrid::smpi::Comm::f2c(*comm)->size();
304   MPI_Datatype* sendtypes = new MPI_Datatype[size];
305   MPI_Datatype* recvtypes = new MPI_Datatype[size];
306   for(int i=0; i< size; i++){
307     if(FORT_IN_PLACE(sendbuf)!=MPI_IN_PLACE)
308       sendtypes[i] = simgrid::smpi::Datatype::f2c(old_sendtypes[i]);
309     recvtypes[i] = simgrid::smpi::Datatype::f2c(old_recvtypes[i]);
310   }
311   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
312  *ierr = MPI_Ialltoallw( sendbuf, sendcnts, sdispls, sendtypes, recvbuf, recvcnts, rdispls,
313                         recvtypes, simgrid::smpi::Comm::f2c(*comm), &req);
314   if(*ierr == MPI_SUCCESS) {
315     *request = req->add_f();
316   }
317   delete[] sendtypes;
318   delete[] recvtypes;
319 }
320
321 void mpi_iexscan_ (void *sendbuf, void *recvbuf, int* count, int* datatype, int* op, int* comm, int* request, int* ierr){
322   MPI_Request req;
323   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
324  *ierr = MPI_Iexscan(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm), &req);
325   if(*ierr == MPI_SUCCESS) {
326     *request = req->add_f();
327   }
328 }
329
330 }