MeDiPack  1.3.0
A Message Differentiation Package
SciComp TU Kaiserslautern
Loading...
Searching...
No Matches
primalFunctions.hpp
Go to the documentation of this file.
1/*
2 * MeDiPack, a Message Differentiation Package
3 *
4 * Copyright (C) 2015-2024 Chair for Scientific Computing (SciComp), University of Kaiserslautern-Landau
5 * Homepage: http://scicomp.rptu.de
6 * Contact: Prof. Nicolas R. Gauger (codi@scicomp.uni-kl.de)
7 *
8 * Lead developers: Max Sagebaum (SciComp, University of Kaiserslautern-Landau)
9 *
10 * This file is part of MeDiPack (http://scicomp.rptu.de/software/codi).
11 *
12 * MeDiPack is free software: you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public
14 * License as published by the Free Software Foundation, either
15 * version 3 of the License, or (at your option) any later version.
16 *
17 * MeDiPack is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
20 *
21 * See the GNU Lesser General Public License for more details.
22 * You should have received a copy of the GNU
23 * Lesser General Public License along with MeDiPack.
24 * If not, see <http://www.gnu.org/licenses/>.
25 *
26 * Authors: Max Sagebaum, Tim Albring (SciComp, University of Kaiserslautern-Landau)
27 */
28
29#pragma once
30
31#include <iostream>
32
33#include "ampiMisc.h"
34#include "async.hpp"
35#include "message.hpp"
36#include "enums.hpp"
38
42namespace medi {
43
44#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
45 template<typename DATATYPE>
46 void AMPI_Send_pri(typename DATATYPE::PrimalType* bufAdjoints, int bufSize, int count, DATATYPE* datatype, int dest, int tag, AMPI_Comm comm) {
47 MEDI_UNUSED(count);
48 MPI_Send(bufAdjoints, bufSize, datatype->getADTool().getPrimalMpiType(), dest, tag, comm);
49 }
50#endif
51
52#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
53 template<typename DATATYPE>
54 void AMPI_Isend_pri(typename DATATYPE::PrimalType* bufAdjoints, int bufSize, int count, DATATYPE* datatype, int dest, int tag, AMPI_Comm comm, AMPI_Request* request) {
55 MEDI_UNUSED(count);
56 MPI_Isend(bufAdjoints, bufSize, datatype->getADTool().getPrimalMpiType(), dest, tag, comm, &request->request);
57 }
58#endif
59
60#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
61 template<typename DATATYPE>
62 void AMPI_Bsend_pri(typename DATATYPE::PrimalType* bufAdjoints, int bufSize, int count, DATATYPE* datatype, int dest, int tag, AMPI_Comm comm) {
63 MEDI_UNUSED(count);
64 MPI_Bsend(bufAdjoints, bufSize, datatype->getADTool().getPrimalMpiType(), dest, tag, comm);
65 }
66#endif
67
68#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
69 template<typename DATATYPE>
70 void AMPI_Ibsend_pri(typename DATATYPE::PrimalType* bufAdjoints, int bufSize, int count, DATATYPE* datatype, int dest, int tag, AMPI_Comm comm, AMPI_Request* request) {
71 MEDI_UNUSED(count);
72 MPI_Ibsend(bufAdjoints, bufSize, datatype->getADTool().getPrimalMpiType(), dest, tag, comm, &request->request);
73 }
74#endif
75
76#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
77 template<typename DATATYPE>
78 void AMPI_Ssend_pri(typename DATATYPE::PrimalType* bufAdjoints, int bufSize, int count, DATATYPE* datatype, int dest, int tag, AMPI_Comm comm) {
79 MEDI_UNUSED(count);
80 MPI_Ssend(bufAdjoints, bufSize, datatype->getADTool().getPrimalMpiType(), dest, tag, comm);
81 }
82#endif
83
84#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
85 template<typename DATATYPE>
86 void AMPI_Issend_pri(typename DATATYPE::PrimalType* bufAdjoints, int bufSize, int count, DATATYPE* datatype, int dest, int tag, AMPI_Comm comm, AMPI_Request* request) {
87 MEDI_UNUSED(count);
88 MPI_Issend(bufAdjoints, bufSize, datatype->getADTool().getPrimalMpiType(), dest, tag, comm, &request->request);
89 }
90#endif
91
92#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
93 template<typename DATATYPE>
94 void AMPI_Rsend_pri(typename DATATYPE::PrimalType* bufAdjoints, int bufSize, int count, DATATYPE* datatype, int dest, int tag, AMPI_Comm comm) {
95 MEDI_UNUSED(count);
96 MPI_Rsend(bufAdjoints, bufSize, datatype->getADTool().getPrimalMpiType(), dest, tag, comm);
97 }
98#endif
99
100#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
101 template<typename DATATYPE>
102 void AMPI_Irsend_pri(typename DATATYPE::PrimalType* bufAdjoints, int bufSize, int count, DATATYPE* datatype, int dest, int tag, AMPI_Comm comm, AMPI_Request* request) {
103 MEDI_UNUSED(count);
104 MPI_Irsend(bufAdjoints, bufSize, datatype->getADTool().getPrimalMpiType(), dest, tag, comm, &request->request);
105 }
106#endif
107
108#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
109 template<typename DATATYPE>
110 void AMPI_Recv_pri(typename DATATYPE::PrimalType* bufAdjoints, int bufSize, int count, DATATYPE* datatype, int src, int tag, AMPI_Comm comm, AMPI_Status* status, RecvAdjCall reverse_call) {
111 MEDI_UNUSED(count);
112 MEDI_UNUSED(reverse_call);
113 MPI_Recv(bufAdjoints, bufSize, datatype->getADTool().getPrimalMpiType(), src, tag, comm, status);
114 }
115#endif
116
117#if MEDI_MPI_VERSION_3_0 <= MEDI_MPI_TARGET
118 template<typename DATATYPE>
119 void AMPI_Mrecv_pri(typename DATATYPE::PrimalType* bufAdjoints, int bufSize, int count, DATATYPE* datatype, AMPI_Message* message, AMPI_Status* status, RecvAdjCall reverse_call) {
120 MEDI_UNUSED(count);
121 MEDI_UNUSED(reverse_call);
122 MPI_Mrecv(bufAdjoints, bufSize, datatype->getADTool().getPrimalMpiType(), &message->message, status);
123 }
124#endif
125
126#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
127 template<typename DATATYPE>
128 void AMPI_Irecv_pri(typename DATATYPE::PrimalType* bufAdjoints, int bufSize, int count, DATATYPE* datatype, int src, int tag, AMPI_Comm comm, AMPI_Request* request, IrecvAdjCall reverse_call) {
129 MEDI_UNUSED(count);
130 MEDI_UNUSED(reverse_call);
131 MPI_Irecv(bufAdjoints, bufSize, datatype->getADTool().getPrimalMpiType(), src, tag, comm, &request->request);
132 }
133#endif
134
135#if MEDI_MPI_VERSION_3_0 <= MEDI_MPI_TARGET
136 template<typename DATATYPE>
137 void AMPI_Imrecv_pri(typename DATATYPE::PrimalType* bufAdjoints, int bufSize, int count, DATATYPE* datatype, AMPI_Message* message, AMPI_Request* request, IrecvAdjCall reverse_call) {
138 MEDI_UNUSED(count);
139 MEDI_UNUSED(reverse_call);
140 MPI_Imrecv(bufAdjoints, bufSize, datatype->getADTool().getPrimalMpiType(), &message->message, &request->request);
141 }
142#endif
143
144#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
145 template<typename SENDTYPE, typename RECVTYPE>
146 void AMPI_Sendrecv_pri(typename SENDTYPE::PrimalType* sendbuf, int sendbufSize, int sendcount, SENDTYPE* sendtype, int dest, int sendtag,
147 typename RECVTYPE::PrimalType* recvbuf, int recvbufSize, int recvcount, RECVTYPE* recvtype, int source, int recvtag, AMPI_Comm comm, AMPI_Status* status) {
148#endif
149
150 MEDI_UNUSED(sendcount);
151 MEDI_UNUSED(recvcount);
152 MPI_Sendrecv(sendbuf, sendbufSize, sendtype->getADTool().getPrimalMpiType(), dest, sendtag, recvbuf, recvbufSize, recvtype->getADTool().getPrimalMpiType(), source, recvtag, comm, status);
153 }
154
155#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
156 template<typename DATATYPE>
157 void AMPI_Bcast_wrap_pri(typename DATATYPE::PrimalType* &sendbufAdjoints, int sendbufSize, typename DATATYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int count, DATATYPE* datatype, int root, AMPI_Comm comm) {
158 MEDI_UNUSED(count);
159
160 if(root == getCommRank(comm)) {
161 std::swap(sendbufAdjoints, recvbufAdjoints);
162 }
163 MPI_Bcast(recvbufAdjoints, recvbufSize, datatype->getADTool().getPrimalMpiType(), root, comm);
164 }
165#endif
166
167#if MEDI_MPI_VERSION_3_0 <= MEDI_MPI_TARGET
168 template<typename DATATYPE>
169 void AMPI_Ibcast_wrap_pri(typename DATATYPE::PrimalType* &sendbufAdjoints, int sendbufSize, typename DATATYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int count, DATATYPE* datatype, int root, AMPI_Comm comm, AMPI_Request* request) {
170 MEDI_UNUSED(count);
171
172 if(root == getCommRank(comm)) {
173 std::swap(sendbufAdjoints, recvbufAdjoints);
174 }
175 MPI_Ibcast(recvbufAdjoints, recvbufSize, datatype->getADTool().getPrimalMpiType(), root, comm, &request->request);
176 }
177#endif
178
179#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
180 template<typename SENDTYPE, typename RECVTYPE>
181 void AMPI_Scatter_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE* recvtype, int root, AMPI_Comm comm) {
182 MEDI_UNUSED(sendcount);
183 MEDI_UNUSED(recvcount);
184
185 MPI_Scatter(sendbufAdjoints, sendbufSize, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufSize, recvtype->getADTool().getPrimalMpiType(), root, comm);
186 }
187#endif
188
189#if MEDI_MPI_VERSION_3_0 <= MEDI_MPI_TARGET
190 template<typename SENDTYPE, typename RECVTYPE>
191 void AMPI_Iscatter_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE* recvtype, int root, AMPI_Comm comm, AMPI_Request* request) {
192 MEDI_UNUSED(sendcount);
193 MEDI_UNUSED(recvcount);
194
195 MPI_Iscatter(sendbufAdjoints, sendbufSize, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufSize, recvtype->getADTool().getPrimalMpiType(), root, comm, &request->request);
196 }
197#endif
198
199#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
200 template<typename SENDTYPE, typename RECVTYPE>
201 void AMPI_Scatterv_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int* sendbufCounts, MEDI_OPTIONAL_CONST int* sendbufDispl, MEDI_OPTIONAL_CONST int* sendcount, MEDI_OPTIONAL_CONST int* displs, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE* recvtype, int root, AMPI_Comm comm) {
202 MEDI_UNUSED(sendcount);
203 MEDI_UNUSED(recvcount);
204
205 MPI_Scatterv(sendbufAdjoints, sendbufCounts, sendbufDispl, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufSize, recvtype->getADTool().getPrimalMpiType(), root, comm);
206 }
207#endif
208
209#if MEDI_MPI_VERSION_3_0 <= MEDI_MPI_TARGET
210 template<typename SENDTYPE, typename RECVTYPE>
211 void AMPI_Iscatterv_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int* sendbufCounts, MEDI_OPTIONAL_CONST int* sendbufDispl, MEDI_OPTIONAL_CONST int* sendcount, MEDI_OPTIONAL_CONST int* displs, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE* recvtype, int root, AMPI_Comm comm, AMPI_Request* request) {
212 MEDI_UNUSED(sendcount);
213 MEDI_UNUSED(recvcount);
214
215 MPI_Iscatterv(sendbufAdjoints, sendbufCounts, sendbufDispl, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufSize, recvtype->getADTool().getPrimalMpiType(), root, comm, &request->request);
216 }
217#endif
218
219#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
220 template<typename SENDTYPE, typename RECVTYPE>
221 void AMPI_Gather_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE* recvtype, int root, AMPI_Comm comm) {
222 MEDI_UNUSED(sendcount);
223 MEDI_UNUSED(recvcount);
224
225 MPI_Gather(sendbufAdjoints, sendbufSize, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufSize, recvtype->getADTool().getPrimalMpiType(), root, comm);
226 }
227#endif
228
229#if MEDI_MPI_VERSION_3_0 <= MEDI_MPI_TARGET
230 template<typename SENDTYPE, typename RECVTYPE>
231 void AMPI_Igather_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE* recvtype, int root, AMPI_Comm comm, AMPI_Request* request) {
232 MEDI_UNUSED(sendcount);
233 MEDI_UNUSED(recvcount);
234
235 MPI_Igather(sendbufAdjoints, sendbufSize, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufSize, recvtype->getADTool().getPrimalMpiType(), root, comm, &request->request);
236 }
237#endif
238
239#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
240 template<typename SENDTYPE, typename RECVTYPE>
241 void AMPI_Gatherv_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int* recvbufCounts, MEDI_OPTIONAL_CONST int* recvbufDispls, MEDI_OPTIONAL_CONST int* recvcounts, MEDI_OPTIONAL_CONST int* displs, RECVTYPE* recvtype, int root, AMPI_Comm comm) {
242 MEDI_UNUSED(sendcount);
243 MEDI_UNUSED(recvcounts);
244 MEDI_UNUSED(displs);
245
246 MPI_Gatherv(sendbufAdjoints, sendbufSize, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufCounts, recvbufDispls, recvtype->getADTool().getPrimalMpiType(), root, comm);
247 }
248#endif
249
250#if MEDI_MPI_VERSION_3_0 <= MEDI_MPI_TARGET
251 template<typename SENDTYPE, typename RECVTYPE>
252 void AMPI_Igatherv_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int* recvbufCounts, MEDI_OPTIONAL_CONST int* recvbufDispls, MEDI_OPTIONAL_CONST int* recvcounts, MEDI_OPTIONAL_CONST int* displs, RECVTYPE* recvtype, int root, AMPI_Comm comm, AMPI_Request* request) {
253 MEDI_UNUSED(sendcount);
254 MEDI_UNUSED(recvcounts);
255 MEDI_UNUSED(displs);
256
257 MPI_Igatherv(sendbufAdjoints, sendbufSize, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufCounts, recvbufDispls, recvtype->getADTool().getPrimalMpiType(), root, comm, &request->request);
258 }
259#endif
260
261#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
262 template<typename SENDTYPE, typename RECVTYPE>
263 void AMPI_Allgather_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE* recvtype, AMPI_Comm comm) {
264 MEDI_UNUSED(sendcount);
265 MEDI_UNUSED(recvcount);
266
267 MPI_Allgather(sendbufAdjoints, sendbufSize, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufSize, recvtype->getADTool().getPrimalMpiType(), comm);
268 }
269#endif
270
271#if MEDI_MPI_VERSION_3_0 <= MEDI_MPI_TARGET
272 template<typename SENDTYPE, typename RECVTYPE>
273 void AMPI_Iallgather_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE* recvtype, AMPI_Comm comm, AMPI_Request* request ) {
274 MEDI_UNUSED(sendcount);
275 MEDI_UNUSED(recvcount);
276
277 MPI_Iallgather(sendbufAdjoints, sendbufSize, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufSize, recvtype->getADTool().getPrimalMpiType(), comm, &request->request);
278 }
279#endif
280
281#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
282 template<typename SENDTYPE, typename RECVTYPE>
283 void AMPI_Allgatherv_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int* recvbufCounts, MEDI_OPTIONAL_CONST int* recvbufDispls, MEDI_OPTIONAL_CONST int* recvcounts, MEDI_OPTIONAL_CONST int* displs, RECVTYPE* recvtype, AMPI_Comm comm) {
284 MEDI_UNUSED(sendcount);
285 MEDI_UNUSED(recvcounts);
286 MEDI_UNUSED(displs);
287
288 MPI_Allgatherv(sendbufAdjoints, sendbufSize, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufCounts, recvbufDispls, recvtype->getADTool().getPrimalMpiType(), comm);
289 }
290#endif
291
292#if MEDI_MPI_VERSION_3_0 <= MEDI_MPI_TARGET
293 template<typename SENDTYPE, typename RECVTYPE>
294 void AMPI_Iallgatherv_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int* recvbufCounts, MEDI_OPTIONAL_CONST int* recvbufDispls, MEDI_OPTIONAL_CONST int* recvcounts, MEDI_OPTIONAL_CONST int* displs, RECVTYPE* recvtype, AMPI_Comm comm, AMPI_Request* request) {
295 MEDI_UNUSED(sendcount);
296 MEDI_UNUSED(recvcounts);
297 MEDI_UNUSED(displs);
298
299 MPI_Iallgatherv(sendbufAdjoints, sendbufSize, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufCounts, recvbufDispls, recvtype->getADTool().getPrimalMpiType(), comm, &request->request);
300 }
301#endif
302
303#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
304 template<typename SENDTYPE, typename RECVTYPE>
305 void AMPI_Alltoall_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE* recvtype, AMPI_Comm comm) {
306 MEDI_UNUSED(sendcount);
307 MEDI_UNUSED(recvcount);
308
309 MPI_Alltoall(sendbufAdjoints, sendbufSize, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufSize, recvtype->getADTool().getPrimalMpiType(), comm);
310 }
311#endif
312
313#if MEDI_MPI_VERSION_3_0 <= MEDI_MPI_TARGET
314 template<typename SENDTYPE, typename RECVTYPE>
315 void AMPI_Ialltoall_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE* recvtype, AMPI_Comm comm, AMPI_Request* request) {
316 MEDI_UNUSED(sendcount);
317 MEDI_UNUSED(recvcount);
318
319 MPI_Ialltoall(sendbufAdjoints, sendbufSize, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufSize, recvtype->getADTool().getPrimalMpiType(), comm, &request->request);
320 }
321#endif
322
323#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
324 template<typename SENDTYPE, typename RECVTYPE>
325 void AMPI_Alltoallv_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int* sendbufCounts, MEDI_OPTIONAL_CONST int* sendbufDispls, MEDI_OPTIONAL_CONST int* sendcounts, MEDI_OPTIONAL_CONST int* sdispls, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int* recvbufCounts, MEDI_OPTIONAL_CONST int* recvbufDispls, MEDI_OPTIONAL_CONST int* recvcounts, MEDI_OPTIONAL_CONST int* rdispls, RECVTYPE* recvtype, AMPI_Comm comm) {
326 MEDI_UNUSED(sendcounts);
327 MEDI_UNUSED(sdispls);
328 MEDI_UNUSED(recvcounts);
329 MEDI_UNUSED(rdispls);
330
331 MPI_Alltoallv(sendbufAdjoints, sendbufCounts, sendbufDispls, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufCounts, recvbufDispls, recvtype->getADTool().getPrimalMpiType(), comm);
332 }
333#endif
334
335#if MEDI_MPI_VERSION_3_0 <= MEDI_MPI_TARGET
336 template<typename SENDTYPE, typename RECVTYPE>
337 void AMPI_Ialltoallv_pri(typename SENDTYPE::PrimalType* &sendbufAdjoints, int* sendbufCounts, MEDI_OPTIONAL_CONST int* sendbufDispls, MEDI_OPTIONAL_CONST int* sendcounts, MEDI_OPTIONAL_CONST int* sdispls, SENDTYPE* sendtype, typename RECVTYPE::PrimalType* &recvbufAdjoints, int* recvbufCounts, MEDI_OPTIONAL_CONST int* recvbufDispls, MEDI_OPTIONAL_CONST int* recvcounts, MEDI_OPTIONAL_CONST int* rdispls, RECVTYPE* recvtype, AMPI_Comm comm, AMPI_Request* request) {
338 MEDI_UNUSED(sendcounts);
339 MEDI_UNUSED(sdispls);
340 MEDI_UNUSED(recvcounts);
341 MEDI_UNUSED(rdispls);
342
343 MPI_Ialltoallv(sendbufAdjoints, sendbufCounts, sendbufDispls, sendtype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufCounts, recvbufDispls, recvtype->getADTool().getPrimalMpiType(), comm, &request->request);
344 }
345#endif
346
347#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
348 template<typename DATATYPE>
349 void AMPI_Reduce_global_pri(typename DATATYPE::PrimalType* &sendbufAdjoints, int sendbufSize, typename DATATYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int count, DATATYPE* datatype, AMPI_Op op, int root, AMPI_Comm comm) {
350 MEDI_UNUSED(count);
351 MEDI_UNUSED(sendbufSize);
352 //TODO: MPI_Reduce(sendbufAdjoints, recvbufAdjoints, recvbufSize, datatype->getADTool().getPrimalMpiType(), TODO, root, comm);
353 std::cout << "Forward reduce not supported." << std::endl;
354 }
355#endif
356
357#if MEDI_MPI_VERSION_3_0 <= MEDI_MPI_TARGET
358 template<typename DATATYPE>
359 void AMPI_Ireduce_global_pri(typename DATATYPE::PrimalType* &sendbufAdjoints, int sendbufSize, typename DATATYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int count, DATATYPE* datatype, AMPI_Op op, int root, AMPI_Comm comm, AMPI_Request* request) {
360 MEDI_UNUSED(count);
361 MEDI_UNUSED(sendbufSize);
362 //TODO: MPI_Ireduce(sendbufAdjoints, recvbufAdjoints, recvbufSize, datatype->getADTool().getPrimalMpiType(), TODO, root, comm, &request->request);
363 std::cout << "Forward reduce not supported." << std::endl;
364 }
365#endif
366
367#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
368 template<typename DATATYPE>
369 void AMPI_Allreduce_global_pri(typename DATATYPE::PrimalType* &sendbufAdjoints, int sendbufSize, typename DATATYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int count, DATATYPE* datatype, AMPI_Op op, AMPI_Comm comm) {
370 MEDI_UNUSED(count);
371 MEDI_UNUSED(sendbufSize);
372
373 //TODO MPI_Allreduce(sendbufAdjoints, datatype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufSize, datatype->getADTool().getPrimalMpiType(), TODO, comm);
374 std::cout << "Forward reduce not supported." << std::endl;
375 }
376#endif
377
378#if MEDI_MPI_VERSION_3_0 <= MEDI_MPI_TARGET
379 template<typename DATATYPE>
380 void AMPI_Iallreduce_global_pri(typename DATATYPE::PrimalType* &sendbufAdjoints, int sendbufSize, typename DATATYPE::PrimalType* &recvbufAdjoints, int recvbufSize, int count, DATATYPE* datatype, AMPI_Op op, AMPI_Comm comm, AMPI_Request* request) {
381 MEDI_UNUSED(op);
382 MEDI_UNUSED(count);
383
384 //TODO MPI_Iallreduce(sendbufAdjoints, datatype->getADTool().getPrimalMpiType(), recvbufAdjoints, recvbufSize, datatype->getADTool().getPrimalMpiType(), comm, &request->request);
385 std::cout << "Forward reduce not supported." << std::endl;
386 }
387#endif
388}
#define AMPI_Status
Definition ampiDefinitions.h:772
#define AMPI_Comm
Definition ampiDefinitions.h:778
#define MEDI_OPTIONAL_CONST
Definition macros.h:70
#define MEDI_UNUSED(name)
Definition macros.h:108
Global namespace for MeDiPack - Message Differentiation Package.
Definition adjointInterface.hpp:37
void AMPI_Ibcast_wrap_pri(typename DATATYPE::PrimalType *&sendbufAdjoints, int sendbufSize, typename DATATYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int count, DATATYPE *datatype, int root, MPI_Comm comm, AMPI_Request *request)
Definition primalFunctions.hpp:169
void AMPI_Isend_pri(typename DATATYPE::PrimalType *bufAdjoints, int bufSize, int count, DATATYPE *datatype, int dest, int tag, MPI_Comm comm, AMPI_Request *request)
Definition primalFunctions.hpp:54
void AMPI_Bcast_wrap_pri(typename DATATYPE::PrimalType *&sendbufAdjoints, int sendbufSize, typename DATATYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int count, DATATYPE *datatype, int root, MPI_Comm comm)
Definition primalFunctions.hpp:157
void AMPI_Iallgather_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE *recvtype, MPI_Comm comm, AMPI_Request *request)
Definition primalFunctions.hpp:273
void AMPI_Issend_pri(typename DATATYPE::PrimalType *bufAdjoints, int bufSize, int count, DATATYPE *datatype, int dest, int tag, MPI_Comm comm, AMPI_Request *request)
Definition primalFunctions.hpp:86
void AMPI_Ireduce_global_pri(typename DATATYPE::PrimalType *&sendbufAdjoints, int sendbufSize, typename DATATYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int count, DATATYPE *datatype, AMPI_Op op, int root, MPI_Comm comm, AMPI_Request *request)
Definition primalFunctions.hpp:359
void AMPI_Irsend_pri(typename DATATYPE::PrimalType *bufAdjoints, int bufSize, int count, DATATYPE *datatype, int dest, int tag, MPI_Comm comm, AMPI_Request *request)
Definition primalFunctions.hpp:102
void AMPI_Rsend_pri(typename DATATYPE::PrimalType *bufAdjoints, int bufSize, int count, DATATYPE *datatype, int dest, int tag, MPI_Comm comm)
Definition primalFunctions.hpp:94
void AMPI_Igather_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE *recvtype, int root, MPI_Comm comm, AMPI_Request *request)
Definition primalFunctions.hpp:231
void AMPI_Bsend_pri(typename DATATYPE::PrimalType *bufAdjoints, int bufSize, int count, DATATYPE *datatype, int dest, int tag, MPI_Comm comm)
Definition primalFunctions.hpp:62
void AMPI_Sendrecv_pri(typename SENDTYPE::PrimalType *sendbuf, int sendbufSize, int sendcount, SENDTYPE *sendtype, int dest, int sendtag, typename RECVTYPE::PrimalType *recvbuf, int recvbufSize, int recvcount, RECVTYPE *recvtype, int source, int recvtag, MPI_Comm comm, MPI_Status *status)
Definition primalFunctions.hpp:146
IrecvAdjCall
Definition enums.hpp:40
void AMPI_Ibsend_pri(typename DATATYPE::PrimalType *bufAdjoints, int bufSize, int count, DATATYPE *datatype, int dest, int tag, MPI_Comm comm, AMPI_Request *request)
Definition primalFunctions.hpp:70
RecvAdjCall
Definition enums.hpp:44
void AMPI_Iallreduce_global_pri(typename DATATYPE::PrimalType *&sendbufAdjoints, int sendbufSize, typename DATATYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int count, DATATYPE *datatype, AMPI_Op op, MPI_Comm comm, AMPI_Request *request)
Definition primalFunctions.hpp:380
void AMPI_Iallgatherv_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int *recvbufCounts, const int *recvbufDispls, const int *recvcounts, const int *displs, RECVTYPE *recvtype, MPI_Comm comm, AMPI_Request *request)
Definition primalFunctions.hpp:294
void AMPI_Alltoallv_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int *sendbufCounts, const int *sendbufDispls, const int *sendcounts, const int *sdispls, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int *recvbufCounts, const int *recvbufDispls, const int *recvcounts, const int *rdispls, RECVTYPE *recvtype, MPI_Comm comm)
Definition primalFunctions.hpp:325
void AMPI_Imrecv_pri(typename DATATYPE::PrimalType *bufAdjoints, int bufSize, int count, DATATYPE *datatype, AMPI_Message *message, AMPI_Request *request, IrecvAdjCall reverse_call)
Definition primalFunctions.hpp:137
void AMPI_Irecv_pri(typename DATATYPE::PrimalType *bufAdjoints, int bufSize, int count, DATATYPE *datatype, int src, int tag, MPI_Comm comm, AMPI_Request *request, IrecvAdjCall reverse_call)
Definition primalFunctions.hpp:128
void AMPI_Allreduce_global_pri(typename DATATYPE::PrimalType *&sendbufAdjoints, int sendbufSize, typename DATATYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int count, DATATYPE *datatype, AMPI_Op op, MPI_Comm comm)
Definition primalFunctions.hpp:369
void AMPI_Alltoall_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE *recvtype, MPI_Comm comm)
Definition primalFunctions.hpp:305
void AMPI_Iscatterv_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int *sendbufCounts, const int *sendbufDispl, const int *sendcount, const int *displs, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE *recvtype, int root, MPI_Comm comm, AMPI_Request *request)
Definition primalFunctions.hpp:211
void AMPI_Scatter_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE *recvtype, int root, MPI_Comm comm)
Definition primalFunctions.hpp:181
void AMPI_Mrecv_pri(typename DATATYPE::PrimalType *bufAdjoints, int bufSize, int count, DATATYPE *datatype, AMPI_Message *message, MPI_Status *status, RecvAdjCall reverse_call)
Definition primalFunctions.hpp:119
void AMPI_Allgather_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE *recvtype, MPI_Comm comm)
Definition primalFunctions.hpp:263
int getCommRank(MPI_Comm comm)
Helper function that gets the own rank number from the communicator.
Definition mpiTools.h:52
void AMPI_Iscatter_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE *recvtype, int root, MPI_Comm comm, AMPI_Request *request)
Definition primalFunctions.hpp:191
void AMPI_Scatterv_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int *sendbufCounts, const int *sendbufDispl, const int *sendcount, const int *displs, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE *recvtype, int root, MPI_Comm comm)
Definition primalFunctions.hpp:201
void AMPI_Recv_pri(typename DATATYPE::PrimalType *bufAdjoints, int bufSize, int count, DATATYPE *datatype, int src, int tag, MPI_Comm comm, MPI_Status *status, RecvAdjCall reverse_call)
Definition primalFunctions.hpp:110
void AMPI_Gatherv_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int *recvbufCounts, const int *recvbufDispls, const int *recvcounts, const int *displs, RECVTYPE *recvtype, int root, MPI_Comm comm)
Definition primalFunctions.hpp:241
void AMPI_Igatherv_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int *recvbufCounts, const int *recvbufDispls, const int *recvcounts, const int *displs, RECVTYPE *recvtype, int root, MPI_Comm comm, AMPI_Request *request)
Definition primalFunctions.hpp:252
void AMPI_Ssend_pri(typename DATATYPE::PrimalType *bufAdjoints, int bufSize, int count, DATATYPE *datatype, int dest, int tag, MPI_Comm comm)
Definition primalFunctions.hpp:78
void AMPI_Allgatherv_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int *recvbufCounts, const int *recvbufDispls, const int *recvcounts, const int *displs, RECVTYPE *recvtype, MPI_Comm comm)
Definition primalFunctions.hpp:283
void AMPI_Send_pri(typename DATATYPE::PrimalType *bufAdjoints, int bufSize, int count, DATATYPE *datatype, int dest, int tag, MPI_Comm comm)
Definition primalFunctions.hpp:46
void AMPI_Gather_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE *recvtype, int root, MPI_Comm comm)
Definition primalFunctions.hpp:221
void AMPI_Ialltoallv_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int *sendbufCounts, const int *sendbufDispls, const int *sendcounts, const int *sdispls, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int *recvbufCounts, const int *recvbufDispls, const int *recvcounts, const int *rdispls, RECVTYPE *recvtype, MPI_Comm comm, AMPI_Request *request)
Definition primalFunctions.hpp:337
void AMPI_Ialltoall_pri(typename SENDTYPE::PrimalType *&sendbufAdjoints, int sendbufSize, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int recvcount, RECVTYPE *recvtype, MPI_Comm comm, AMPI_Request *request)
Definition primalFunctions.hpp:315
void AMPI_Reduce_global_pri(typename DATATYPE::PrimalType *&sendbufAdjoints, int sendbufSize, typename DATATYPE::PrimalType *&recvbufAdjoints, int recvbufSize, int count, DATATYPE *datatype, AMPI_Op op, int root, MPI_Comm comm)
Definition primalFunctions.hpp:349
Stores additional information for a MPI_Message.
Definition message.hpp:44
MPI_Message message
Definition message.hpp:45
Structure for the special handling of the MPI_Op structure.
Definition op.hpp:50
Definition async.hpp:45
MPI_Request request
Definition async.hpp:46