1 /* Copyright (c) 2010-2017. The SimGrid Team. All rights reserved. */
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. */
7 #include "smpi_coll.hpp"
8 #include "smpi_comm.hpp"
9 #include "smpi_datatype.hpp"
10 #include "smpi_op.hpp"
12 extern "C" { // This should really use the C linkage to be usable from Fortran
14 void mpi_barrier_(int* comm, int* ierr) {
15 *ierr = MPI_Barrier(simgrid::smpi::Comm::f2c(*comm));
18 void mpi_bcast_(void *buf, int* count, int* datatype, int* root, int* comm, int* ierr) {
19 *ierr = MPI_Bcast(buf, *count, simgrid::smpi::Datatype::f2c(*datatype), *root, simgrid::smpi::Comm::f2c(*comm));
22 void mpi_reduce_(void* sendbuf, void* recvbuf, int* count, int* datatype, int* op, int* root, int* comm, int* ierr) {
23 sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
24 sendbuf = static_cast<char *>( FORT_BOTTOM(sendbuf));
25 recvbuf = static_cast<char *>( FORT_BOTTOM(recvbuf));
26 *ierr = MPI_Reduce(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), *root, simgrid::smpi::Comm::f2c(*comm));
29 void mpi_allreduce_(void* sendbuf, void* recvbuf, int* count, int* datatype, int* op, int* comm, int* ierr) {
30 sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
31 *ierr = MPI_Allreduce(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm));
34 void mpi_reduce_scatter_(void* sendbuf, void* recvbuf, int* recvcounts, int* datatype, int* op, int* comm, int* ierr) {
35 sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
36 *ierr = MPI_Reduce_scatter(sendbuf, recvbuf, recvcounts, simgrid::smpi::Datatype::f2c(*datatype),
37 simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm));
40 void mpi_scatter_(void* sendbuf, int* sendcount, int* sendtype, void* recvbuf, int* recvcount, int* recvtype,
41 int* root, int* comm, int* ierr) {
42 recvbuf = static_cast<char *>( FORT_IN_PLACE(recvbuf));
43 *ierr = MPI_Scatter(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
44 recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm));
47 void mpi_scatterv_(void* sendbuf, int* sendcounts, int* displs, int* sendtype,
48 void* recvbuf, int* recvcount, int* recvtype, int* root, int* comm, int* ierr) {
49 recvbuf = static_cast<char *>( FORT_IN_PLACE(recvbuf));
50 *ierr = MPI_Scatterv(sendbuf, sendcounts, displs, simgrid::smpi::Datatype::f2c(*sendtype),
51 recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm));
54 void mpi_gather_(void* sendbuf, int* sendcount, int* sendtype, void* recvbuf, int* recvcount, int* recvtype,
55 int* root, int* comm, int* ierr) {
56 sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
57 sendbuf = sendbuf!=MPI_IN_PLACE ? static_cast<char *>( FORT_BOTTOM(sendbuf)) : MPI_IN_PLACE;
58 recvbuf = static_cast<char *>( FORT_BOTTOM(recvbuf));
59 *ierr = MPI_Gather(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
60 recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm));
63 void mpi_gatherv_(void* sendbuf, int* sendcount, int* sendtype,
64 void* recvbuf, int* recvcounts, int* displs, int* recvtype, int* root, int* comm, int* ierr) {
65 sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
66 sendbuf = sendbuf!=MPI_IN_PLACE ? static_cast<char *>( FORT_BOTTOM(sendbuf)) : MPI_IN_PLACE;
67 recvbuf = static_cast<char *>( FORT_BOTTOM(recvbuf));
68 *ierr = MPI_Gatherv(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
69 recvbuf, recvcounts, displs, simgrid::smpi::Datatype::f2c(*recvtype), *root, simgrid::smpi::Comm::f2c(*comm));
72 void mpi_allgather_(void* sendbuf, int* sendcount, int* sendtype, void* recvbuf, int* recvcount, int* recvtype,
73 int* comm, int* ierr) {
74 sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
75 *ierr = MPI_Allgather(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
76 recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm));
79 void mpi_allgatherv_(void* sendbuf, int* sendcount, int* sendtype,
80 void* recvbuf, int* recvcounts,int* displs, int* recvtype, int* comm, int* ierr) {
81 sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
82 *ierr = MPI_Allgatherv(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
83 recvbuf, recvcounts, displs, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm));
86 void mpi_scan_(void* sendbuf, void* recvbuf, int* count, int* datatype, int* op, int* comm, int* ierr) {
87 *ierr = MPI_Scan(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype),
88 simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm));
91 void mpi_alltoall_(void* sendbuf, int* sendcount, int* sendtype,
92 void* recvbuf, int* recvcount, int* recvtype, int* comm, int* ierr) {
93 *ierr = MPI_Alltoall(sendbuf, *sendcount, simgrid::smpi::Datatype::f2c(*sendtype),
94 recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm));
97 void mpi_alltoallv_(void* sendbuf, int* sendcounts, int* senddisps, int* sendtype,
98 void* recvbuf, int* recvcounts, int* recvdisps, int* recvtype, int* comm, int* ierr) {
99 *ierr = MPI_Alltoallv(sendbuf, sendcounts, senddisps, simgrid::smpi::Datatype::f2c(*sendtype),
100 recvbuf, recvcounts, recvdisps, simgrid::smpi::Datatype::f2c(*recvtype), simgrid::smpi::Comm::f2c(*comm));
103 void mpi_reduce_local_ (void *inbuf, void *inoutbuf, int* count, int* datatype, int* op, int* ierr){
105 *ierr = MPI_Reduce_local(inbuf, inoutbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op));
108 void mpi_reduce_scatter_block_ (void *sendbuf, void *recvbuf, int* recvcount, int* datatype, int* op, int* comm,
111 sendbuf = static_cast<char *>( FORT_IN_PLACE(sendbuf));
112 *ierr = MPI_Reduce_scatter_block(sendbuf, recvbuf, *recvcount, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op),
113 simgrid::smpi::Comm::f2c(*comm));
116 void mpi_alltoallw_ ( void *sendbuf, int *sendcnts, int *sdispls, int* sendtypes, void *recvbuf, int *recvcnts,
117 int *rdispls, int* recvtypes, int* comm, int* ierr){
118 *ierr = MPI_Alltoallw( sendbuf, sendcnts, sdispls, reinterpret_cast<MPI_Datatype*>(sendtypes), recvbuf, recvcnts, rdispls,
119 reinterpret_cast<MPI_Datatype*>(recvtypes), simgrid::smpi::Comm::f2c(*comm));
122 void mpi_exscan_ (void *sendbuf, void *recvbuf, int* count, int* datatype, int* op, int* comm, int* ierr){
123 *ierr = MPI_Exscan(sendbuf, recvbuf, *count, simgrid::smpi::Datatype::f2c(*datatype), simgrid::smpi::Op::f2c(*op), simgrid::smpi::Comm::f2c(*comm));