Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
Add new entry in Release_Notes.
[simgrid.git] / src / smpi / bindings / smpi_f77_coll.cpp
1 /* Copyright (c) 2010-2023. 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   std::vector<MPI_Datatype> sendtypes(size);
123   std::vector<MPI_Datatype> recvtypes(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.data(), recvbuf, recvcnts, rdispls, recvtypes.data(),
131                         simgrid::smpi::Comm::f2c(*comm));
132 }
133
134 void mpi_exscan_ (void *sendbuf, void *recvbuf, int* count, int* datatype, int* op, int* comm, int* ierr){
135  *ierr = MPI_Exscan(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm));
136 }
137
138 void mpi_ibarrier_(int* comm, int* request, int* ierr) {
139   MPI_Request req;
140   *ierr = MPI_Ibarrier(simgrid::smpi::Comm::f2c(*comm), &req);
141   if(*ierr == MPI_SUCCESS) {
142     *request = req->c2f();
143   }
144 }
145
146 void mpi_ibcast_(void *buf, int* count, int* datatype, int* root, int* comm, int* request, int* ierr) {
147   MPI_Request req;
148   *ierr = MPI_Ibcast(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *root, simgrid::smpi::Comm::f2c(*comm), &req);
149   if(*ierr == MPI_SUCCESS) {
150     *request = req->c2f();
151   }
152 }
153
154 void mpi_ireduce_(void* sendbuf, void* recvbuf, int* count, int* datatype, int* op, int* root, int* comm, int* request, int* ierr) {
155   MPI_Request req;
156   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
157   sendbuf = static_cast<char *>( FORT_BOTTOM(sendbuf));
158   recvbuf = static_cast<char *>( FORT_BOTTOM(recvbuf));
159   *ierr = MPI_Ireduce(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), *root, simgrid::smpi::Comm::f2c(*comm), &req);
160   if(*ierr == MPI_SUCCESS) {
161     *request = req->c2f();
162   }
163 }
164
165 void mpi_iallreduce_(void* sendbuf, void* recvbuf, int* count, int* datatype, int* op, int* comm, int* request, int* ierr) {
166   MPI_Request req;
167   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
168   *ierr = MPI_Iallreduce(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm), &req);
169   if(*ierr == MPI_SUCCESS) {
170     *request = req->c2f();
171   }
172 }
173
174 void mpi_ireduce_scatter_(void* sendbuf, void* recvbuf, int* recvcounts, int* datatype, int* op, int* comm, int* request, int* ierr) {
175   MPI_Request req;
176   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
177   *ierr = MPI_Ireduce_scatter(sendbuf, recvbuf, recvcounts, simgrid::smpi::Datatype::f2c(*datatype),
178                         simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm), &req);
179   if(*ierr == MPI_SUCCESS) {
180     *request = req->c2f();
181   }
182 }
183
184 void mpi_iscatter_(void* sendbuf, int* sendcount, int* sendtype, void* recvbuf, int* recvcount, int* recvtype,
185                    int* root, int* comm, int* request, int* ierr) {
186   MPI_Request req;
187   recvbuf = static_cast<char *>( FORT_IN_PLACE(recvbuf));
188   *ierr = MPI_Iscatter(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
189                       recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm), &req);
190   if(*ierr == MPI_SUCCESS) {
191     *request = req->c2f();
192   }
193 }
194
195 void mpi_iscatterv_(void* sendbuf, int* sendcounts, int* displs, int* sendtype,
196                    void* recvbuf, int* recvcount, int* recvtype, int* root, int* comm, int* request, int* ierr) {
197   MPI_Request req;
198   recvbuf = static_cast<char *>( FORT_IN_PLACE(recvbuf));
199   *ierr = MPI_Iscatterv(sendbuf, sendcounts, displs, simgrid::smpi::Datatype::f2c(*sendtype),
200                       recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm), &req);
201   if(*ierr == MPI_SUCCESS) {
202     *request = req->c2f();
203   }
204 }
205
206 void mpi_igather_(void* sendbuf, int* sendcount, int* sendtype, void* recvbuf, int* recvcount, int* recvtype,
207                   int* root, int* comm, int* request, int* ierr) {
208   MPI_Request req;
209   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
210   sendbuf = sendbuf!=MPI_IN_PLACE ? static_cast<char *>( FORT_BOTTOM(sendbuf)) : MPI_IN_PLACE;
211   recvbuf = static_cast<char *>( FORT_BOTTOM(recvbuf));
212   *ierr = MPI_Igather(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
213                      recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm), &req);
214   if(*ierr == MPI_SUCCESS) {
215     *request = req->c2f();
216   }
217 }
218
219 void mpi_igatherv_(void* sendbuf, int* sendcount, int* sendtype,
220                   void* recvbuf, int* recvcounts, int* displs, int* recvtype, int* root, int* comm, int* request, int* ierr) {
221   MPI_Request req;
222   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
223   sendbuf = sendbuf!=MPI_IN_PLACE ? static_cast<char *>( FORT_BOTTOM(sendbuf)) : MPI_IN_PLACE;
224   recvbuf = static_cast<char *>( FORT_BOTTOM(recvbuf));
225   *ierr = MPI_Igatherv(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
226                      recvbuf, recvcounts, displs, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm), &req);
227   if(*ierr == MPI_SUCCESS) {
228     *request = req->c2f();
229   }
230 }
231
232 void mpi_iallgather_(void* sendbuf, int* sendcount, int* sendtype, void* recvbuf, int* recvcount, int* recvtype,
233                      int* comm, int* request, int* ierr) {
234   MPI_Request req;
235   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
236   *ierr = MPI_Iallgather(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
237                         recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm), &req);
238   if(*ierr == MPI_SUCCESS) {
239     *request = req->c2f();
240   }
241 }
242
243 void mpi_iallgatherv_(void* sendbuf, int* sendcount, int* sendtype,
244                      void* recvbuf, int* recvcounts,int* displs, int* recvtype, int* comm, int* request, int* ierr) {
245   MPI_Request req;
246   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
247   *ierr = MPI_Iallgatherv(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
248                         recvbuf, recvcounts, displs, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm), &req);
249   if(*ierr == MPI_SUCCESS) {
250     *request = req->c2f();
251   }
252 }
253
254 void mpi_iscan_(void* sendbuf, void* recvbuf, int* count, int* datatype, int* op, int* comm, int* request, int* ierr) {
255   MPI_Request req;
256   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
257   *ierr = MPI_Iscan(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype),
258                    simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm), &req);
259   if(*ierr == MPI_SUCCESS) {
260     *request = req->c2f();
261   }
262 }
263
264 void mpi_ialltoall_(void* sendbuf, int* sendcount, int* sendtype,
265                     void* recvbuf, int* recvcount, int* recvtype, int* comm, int* request, int* ierr) {
266   MPI_Request req;
267   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
268   *ierr = MPI_Ialltoall(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
269                        recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm), &req);
270   if(*ierr == MPI_SUCCESS) {
271     *request = req->c2f();
272   }
273 }
274
275 void mpi_ialltoallv_(void* sendbuf, int* sendcounts, int* senddisps, int* sendtype,
276                     void* recvbuf, int* recvcounts, int* recvdisps, int* recvtype, int* comm, int* request, int* ierr) {
277   MPI_Request req;
278   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
279   *ierr = MPI_Ialltoallv(sendbuf, sendcounts, senddisps, simgrid::smpi::Datatype::f2c(*sendtype),
280                        recvbuf, recvcounts, recvdisps, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm), &req);
281   if(*ierr == MPI_SUCCESS) {
282     *request = req->c2f();
283   }
284 }
285
286 void mpi_ireduce_scatter_block_ (void *sendbuf, void *recvbuf, int* recvcount, int* datatype, int* op, int* comm,
287                                  int* request, int* ierr)
288 {
289   MPI_Request req;
290   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
291  *ierr = MPI_Ireduce_scatter_block(sendbuf, recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op),
292                                   simgrid::smpi::Comm::f2c(*comm), &req);
293   if(*ierr == MPI_SUCCESS) {
294     *request = req->c2f();
295   }
296 }
297
298 void mpi_ialltoallw_ ( void *sendbuf, int *sendcnts, int *sdispls, int* old_sendtypes, void *recvbuf, int *recvcnts,
299                       int *rdispls, int* old_recvtypes, int* comm, int* request, int* ierr){
300   MPI_Request req;
301   int size = simgrid::smpi::Comm::f2c(*comm)->size();
302   std::vector<MPI_Datatype> sendtypes(size);
303   std::vector<MPI_Datatype> recvtypes(size);
304   for(int i=0; i< size; i++){
305     if(FORT_IN_PLACE(sendbuf)!=MPI_IN_PLACE)
306       sendtypes[i] = simgrid::smpi::Datatype::f2c(old_sendtypes[i]);
307     recvtypes[i] = simgrid::smpi::Datatype::f2c(old_recvtypes[i]);
308   }
309   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
310   *ierr   = MPI_Ialltoallw(sendbuf, sendcnts, sdispls, sendtypes.data(), recvbuf, recvcnts, rdispls, recvtypes.data(),
311                          simgrid::smpi::Comm::f2c(*comm), &req);
312   if(*ierr == MPI_SUCCESS) {
313     *request = req->c2f();
314   }
315 }
316
317 void mpi_iexscan_ (void *sendbuf, void *recvbuf, int* count, int* datatype, int* op, int* comm, int* request, int* ierr){
318   MPI_Request req;
319   sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
320  *ierr = MPI_Iexscan(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm), &req);
321   if(*ierr == MPI_SUCCESS) {
322     *request = req->c2f();
323   }
324 }
325
326 }