Logo AND Algorithmique Numérique Distribuée

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