Add const to thread-MPI functions
[gromacs.git] / src / external / thread_mpi / src / p2p_send_recv.c
blobd3a2dcad7e91633690595747c0214b6d830471a3
1 /*
2 This source code file is part of thread_mpi.
3 Written by Sander Pronk, Erik Lindahl, and possibly others.
5 Copyright (c) 2009,2016, Sander Pronk, Erik Lindahl.
6 All rights reserved.
8 Redistribution and use in source and binary forms, with or without
9 modification, are permitted provided that the following conditions are met:
10 1) Redistributions of source code must retain the above copyright
11 notice, this list of conditions and the following disclaimer.
12 2) Redistributions in binary form must reproduce the above copyright
13 notice, this list of conditions and the following disclaimer in the
14 documentation and/or other materials provided with the distribution.
15 3) Neither the name of the copyright holders nor the
16 names of its contributors may be used to endorse or promote products
17 derived from this software without specific prior written permission.
19 THIS SOFTWARE IS PROVIDED BY US ''AS IS'' AND ANY
20 EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
21 WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
22 DISCLAIMED. IN NO EVENT SHALL WE BE LIABLE FOR ANY
23 DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
24 (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
25 LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND
26 ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
27 (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
28 SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
30 If you want to redistribute modifications, please consider that
31 scientific software is very special. Version control is crucial -
32 bugs must be traceable. We will be happy to consider code for
33 inclusion in the official distribution, but derived work should not
34 be called official thread_mpi. Details are found in the README & COPYING
35 files.
38 #ifdef HAVE_TMPI_CONFIG_H
39 #include "tmpi_config.h"
40 #endif
42 #ifdef HAVE_CONFIG_H
43 #include "config.h"
44 #endif
47 #ifdef HAVE_UNISTD_H
48 #include <unistd.h>
49 #endif
51 #include <errno.h>
52 #include <stdlib.h>
53 #include <stdio.h>
54 #include <stdarg.h>
55 #include <string.h>
58 #include "impl.h"
59 #include "p2p.h"
62 /* point-to-point communication exported functions */
64 int tMPI_Send(const void* buf, int count, tMPI_Datatype datatype, int dest,
65 int tag, tMPI_Comm comm)
67 struct envelope *sev;
68 struct tmpi_thread *send_dst;
69 struct tmpi_thread *cur = tMPI_Get_current();
70 struct tmpi_req_ req;
72 #ifdef TMPI_PROFILE
73 tMPI_Profile_count_start(cur);
74 #endif
75 #ifdef TMPI_TRACE
76 tMPI_Trace_print("tMPI_Send(%p, %d, %p, %d, %d, %p)", buf, count,
77 datatype, dest, tag, comm);
78 #endif
79 if (!comm)
81 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_COMM);
83 send_dst = tMPI_Get_thread(comm, dest);
84 if (!send_dst)
86 return tMPI_Error(comm, TMPI_ERR_SEND_DEST);
89 sev = tMPI_Post_send(cur, comm, send_dst, (void*)buf, count, datatype, tag, FALSE);
90 if (sev == NULL)
92 return TMPI_ERR_ENVELOPES;
94 tMPI_Req_init(&req, sev);
95 tMPI_Wait_single(cur, &req);
97 #ifdef TMPI_PROFILE
98 tMPI_Profile_count_stop(cur, TMPIFN_Send);
99 #endif
100 return req.error;
106 int tMPI_Recv(void* buf, int count, tMPI_Datatype datatype, int source,
107 int tag, tMPI_Comm comm, tMPI_Status *status)
109 struct envelope *rev;
110 struct tmpi_thread *recv_src = 0;
111 struct tmpi_thread *cur = tMPI_Get_current();
112 struct tmpi_req_ req;
114 #ifdef TMPI_PROFILE
115 tMPI_Profile_count_start(cur);
116 #endif
117 #ifdef TMPI_TRACE
118 tMPI_Trace_print("tMPI_Recv(%p, %d, %p, %d, %d, %p, %p)", buf, count,
119 datatype, source, tag, comm, status);
120 #endif
121 if (!comm)
123 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_COMM);
126 if (source != TMPI_ANY_SOURCE)
128 recv_src = tMPI_Get_thread(comm, source);
129 if (!recv_src)
131 return tMPI_Error(comm, TMPI_ERR_RECV_SRC);
135 rev = tMPI_Post_match_recv(cur, comm, recv_src, buf, count, datatype, tag,
136 FALSE);
137 if (rev == NULL)
139 return TMPI_ERR_ENVELOPES;
141 tMPI_Req_init(&req, rev);
142 tMPI_Wait_single(cur, &req);
144 tMPI_Set_status(&req, status);
145 #ifdef TMPI_PROFILE
146 tMPI_Profile_count_stop(cur, TMPIFN_Recv);
147 #endif
148 return req.error;
155 int tMPI_Sendrecv(const void *sendbuf, int sendcount, tMPI_Datatype sendtype,
156 int dest, int sendtag, void *recvbuf, int recvcount,
157 tMPI_Datatype recvtype, int source, int recvtag,
158 tMPI_Comm comm, tMPI_Status *status)
160 struct envelope *rev, *sev;
161 struct tmpi_thread *cur = tMPI_Get_current();
162 struct tmpi_thread *recv_src = 0;
163 struct tmpi_thread *send_dst;
164 struct tmpi_req_ sreq, rreq;
165 int ret = TMPI_SUCCESS;
167 #ifdef TMPI_PROFILE
168 tMPI_Profile_count_start(cur);
169 #endif
170 #ifdef TMPI_TRACE
171 tMPI_Trace_print("tMPI_Sendrecv(%p, %d, %p, %d, %d, %p, %d, %p, %d, %d, %p, %p)",
172 sendbuf, sendcount, sendtype, dest, sendtag, recvbuf,
173 recvcount, recvtype, source, recvtag, comm, status);
174 #endif
175 if (!comm)
177 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_COMM);
179 send_dst = tMPI_Get_thread(comm, dest);
180 if (!send_dst)
182 return tMPI_Error(comm, TMPI_ERR_SEND_DEST);
184 if (source != TMPI_ANY_SOURCE)
186 recv_src = tMPI_Get_thread(comm, source);
187 if (!recv_src)
189 return tMPI_Error(comm, TMPI_ERR_RECV_SRC);
193 /* we first prepare to send */
194 sev = tMPI_Post_send(cur, comm, send_dst, (void*)sendbuf, sendcount,
195 sendtype, sendtag, FALSE);
196 if (sev == NULL)
198 return TMPI_ERR_ENVELOPES;
200 tMPI_Req_init(&sreq, sev);
201 /* the we prepare to receive */
202 rev = tMPI_Post_match_recv(cur, comm, recv_src, recvbuf, recvcount,
203 recvtype, recvtag, FALSE);
204 if (rev == NULL)
206 return TMPI_ERR_ENVELOPES;
208 tMPI_Req_init(&rreq, rev);
210 /* fix the pointers */
211 sreq.next = &rreq;
212 sreq.prev = NULL;
213 rreq.prev = &sreq;
214 rreq.next = NULL;
216 /* and wait for our requests */
219 if (tMPI_Test_multi(cur, &sreq, NULL))
221 break;
223 tMPI_Wait_process_incoming(cur);
225 while (TRUE);
227 #ifdef TMPI_PROFILE
228 tMPI_Profile_count_stop(cur, TMPIFN_Sendrecv);
229 #endif
231 tMPI_Set_status(&rreq, status);
232 ret = sreq.error;
233 if (rreq.error != TMPI_SUCCESS)
235 ret = rreq.error;
238 return ret;
245 /* async */
247 int tMPI_Isend(const void* buf, int count, tMPI_Datatype datatype, int dest,
248 int tag, tMPI_Comm comm, tMPI_Request *request)
250 struct tmpi_thread *cur = tMPI_Get_current();
251 struct req_list *rql = &(cur->rql);
252 struct tmpi_req_ *rq = tMPI_Get_req(rql);
253 struct tmpi_thread *send_dst;
254 struct envelope *ev;
256 #ifdef TMPI_PROFILE
257 tMPI_Profile_count_start(cur);
258 #endif
259 #ifdef TMPI_TRACE
260 tMPI_Trace_print("tMPI_Isend(%p, %d, %p, %d, %d, %p, %p)", buf, count,
261 datatype, dest, tag, comm, request);
262 #endif
263 if (!comm)
265 tMPI_Return_req(rql, rq);
266 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_COMM);
268 send_dst = tMPI_Get_thread(comm, dest);
269 if (!send_dst)
271 tMPI_Return_req(rql, rq);
272 return tMPI_Error(comm, TMPI_ERR_SEND_DEST);
274 ev = tMPI_Post_send(cur, comm, send_dst, (void*)buf, count, datatype, tag, TRUE);
275 if (ev == NULL)
277 return TMPI_ERR_ENVELOPES;
279 tMPI_Req_init(rq, ev);
280 *request = rq;
282 #ifdef TMPI_PROFILE
283 tMPI_Profile_count_stop(cur, TMPIFN_Isend);
284 #endif
285 return ev->error;
289 int tMPI_Irecv(void* buf, int count, tMPI_Datatype datatype, int source,
290 int tag, tMPI_Comm comm, tMPI_Request *request)
292 struct tmpi_thread *cur = tMPI_Get_current();
293 struct req_list *rql = &(cur->rql);
294 struct tmpi_req_ *rq = tMPI_Get_req(rql);
295 struct tmpi_thread *recv_src = 0;
296 struct envelope *ev;
298 #ifdef TMPI_PROFILE
299 tMPI_Profile_count_start(cur);
300 #endif
301 #ifdef TMPI_TRACE
302 tMPI_Trace_print("tMPI_Irecv(%p, %d, %p, %d, %d, %p, %p)", buf, count,
303 datatype, source, tag, comm, request);
304 #endif
305 if (!comm)
307 tMPI_Return_req(rql, rq);
308 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_COMM);
311 if (source != TMPI_ANY_SOURCE)
313 recv_src = tMPI_Get_thread(comm, source);
314 if (!recv_src)
316 tMPI_Return_req(rql, rq);
317 return tMPI_Error(comm, TMPI_ERR_RECV_SRC);
320 ev = tMPI_Post_match_recv(cur, comm, recv_src, buf, count, datatype, tag,
321 TRUE);
322 if (ev == NULL)
324 return TMPI_ERR_ENVELOPES;
326 tMPI_Req_init(rq, ev);
327 *request = rq;
328 #ifdef TMPI_PROFILE
329 tMPI_Profile_count_stop(cur, TMPIFN_Irecv);
330 #endif
331 return ev->error;