Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
b21ccbf4f867274314f700a8c85c7a3c19cb8857
[simgrid.git] / src / smpi / mpi.c
1 #include "private.h"
2
3 int MPI_Init(int *argc, char ***argv)
4 {
5         smpi_mpi_init();
6         smpi_bench_begin();
7         return MPI_SUCCESS;
8 }
9
10 int MPI_Finalize()
11 {
12         smpi_bench_end();
13         smpi_mpi_finalize();
14         return MPI_SUCCESS;
15 }
16
17 // right now this just exits the current node, should send abort signal to all
18 // hosts in the communicator;
19 int MPI_Abort(MPI_Comm comm, int errorcode)
20 {
21         smpi_exit(errorcode);
22
23         return 0;
24 }
25
26 int MPI_Comm_size(MPI_Comm comm, int *size)
27 {
28         int retval = MPI_SUCCESS;
29
30         smpi_bench_end();
31
32         if (NULL == comm) {
33                 retval = MPI_ERR_COMM;
34         } else if (NULL == size) {
35                 retval = MPI_ERR_ARG;
36         } else {
37                 *size = comm->size;
38         }
39
40         smpi_bench_begin();
41
42         return retval;
43 }
44
45 int MPI_Comm_rank(MPI_Comm comm, int *rank)
46 {
47         int retval = MPI_SUCCESS;
48
49         smpi_bench_end();
50
51         if (NULL == comm) {
52                 retval = MPI_ERR_COMM;
53         } else if (NULL == rank) {
54                 retval = MPI_ERR_ARG;
55         } else {
56                 *rank = smpi_comm_rank(comm, SIMIX_host_self());
57         }
58
59         smpi_bench_begin();
60
61         return retval;
62 }
63
64 int MPI_Type_size(MPI_Datatype datatype, size_t *size)
65 {
66         int retval = MPI_SUCCESS;
67
68         smpi_bench_end();
69
70         if (NULL == datatype) {
71                 retval = MPI_ERR_TYPE;
72         } else if (NULL == size) {
73                 retval = MPI_ERR_ARG;
74         } else {
75                 *size = datatype->size;
76         }
77
78         smpi_bench_begin();
79
80         return retval;
81 }
82
83 int MPI_Barrier(MPI_Comm comm)
84 {
85         smpi_bench_end();
86         smpi_barrier(comm);
87         smpi_bench_begin();
88         return MPI_SUCCESS;
89 }
90
91 int MPI_Irecv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Request *request)
92 {
93         int retval = MPI_SUCCESS;
94         int dst;
95         smpi_bench_end();
96         dst = smpi_mpi_comm_rank_self(comm);
97         retval = smpi_create_request(buf, count, datatype, src, dst, tag, comm, request);
98         if (NULL != *request) {
99                 smpi_irecv(*request);
100         }
101         smpi_bench_begin();
102         return retval;
103 }
104
105 int MPI_Recv(void *buf, int count, MPI_Datatype datatype, int src, int tag, MPI_Comm comm, MPI_Status *status)
106 {
107         int retval = MPI_SUCCESS;
108         int dst;
109         smpi_mpi_request_t *request;
110         smpi_bench_end();
111         dst = smpi_mpi_comm_rank_self(comm);
112         retval = smpi_create_request(buf, count, datatype, src, dst, tag, comm, &request);
113         if (NULL != request) {
114                 smpi_irecv(request);
115                 smpi_wait(request, status);
116                 // FIXME: mallocator
117                 //xbt_free(request);
118         }
119         smpi_bench_begin();
120         return retval;
121 }
122
123 int MPI_Isend(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm, MPI_Request *request)
124 {
125         int retval = MPI_SUCCESS;
126         int src;
127         smpi_bench_end();
128         src = smpi_mpi_comm_rank_self(comm);
129         retval = smpi_create_request(buf, count, datatype, src, dst, tag, comm, request);
130         if (NULL != *request) {
131                 smpi_isend(*request);
132         }
133         smpi_bench_begin();
134         return retval;
135 }
136
137 int MPI_Send(void *buf, int count, MPI_Datatype datatype, int dst, int tag, MPI_Comm comm)
138 {
139         int retval = MPI_SUCCESS;
140         int src;
141         smpi_mpi_request_t *request;
142         smpi_bench_end();
143         src = smpi_mpi_comm_rank_self(comm);
144         retval = smpi_create_request(buf, count, datatype, src, dst, tag, comm, &request);
145         if (NULL != request) {
146                 smpi_isend(request);
147                 smpi_wait(request, MPI_STATUS_IGNORE);
148                 // FIXME: mallocator
149                 //xbt_free(request)
150         }
151         smpi_bench_begin();
152         return retval;
153 }