Logo AND Algorithmique Numérique Distribuée

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