MeDiPack  1.3.0
A Message Differentiation Package
SciComp TU Kaiserslautern
Loading...
Searching...
No Matches
wrappers.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 "async.hpp"
32#include "ampiMisc.h"
33#include "inPlace.hpp"
34#include "../mpiTools.h"
35
37
41namespace medi {
42
43#if MEDI_MPI_VERSION_1_0 <= MEDI_MPI_TARGET
44 template<typename DATATYPE>
45 int AMPI_Bcast_wrap(typename DATATYPE::Type* bufferSend, typename DATATYPE::Type* bufferRecv, int count, DATATYPE* datatype, int root, AMPI_Comm comm);
46
47 template<typename DATATYPE>
48 inline int AMPI_Bcast(typename DATATYPE::Type* buffer, int count, DATATYPE* datatype, int root, AMPI_Comm comm) {
49 return AMPI_Bcast_wrap<DATATYPE>(AMPI_IN_PLACE, buffer, count, datatype, root, comm);
50 }
51
52 inline int MPI_Bcast_wrap(void* bufferSend, void* bufferRecv, int count, MPI_Datatype type, int root, MPI_Comm comm) {
53 MEDI_UNUSED(bufferSend);
54 return MPI_Bcast(bufferRecv, count, type, root, comm);
55 }
56#endif
57
58#if MEDI_MPI_VERSION_3_0 <= MEDI_MPI_TARGET
59
60 template<typename DATATYPE>
61 int AMPI_Ibcast_wrap(typename DATATYPE::Type* bufferSend, typename DATATYPE::Type* bufferRecv, int count, DATATYPE* datatype, int root, AMPI_Comm comm, AMPI_Request* request);
62
63 template<typename DATATYPE>
64 inline int AMPI_Ibcast(typename DATATYPE::Type* buffer, int count, DATATYPE* datatype, int root, AMPI_Comm comm, AMPI_Request* request) {
65 return AMPI_Ibcast_wrap<DATATYPE>(AMPI_IN_PLACE, buffer, count, datatype, root, comm, request);
66 }
67
68 inline int MPI_Ibcast_wrap(void* bufferSend, void* bufferRecv, int count, MPI_Datatype type, int root, MPI_Comm comm, MPI_Request* request) {
69 MEDI_UNUSED(bufferSend);
70 return MPI_Ibcast(bufferRecv, count, type, root, comm, request);
71 }
72#endif
73
74 template<typename DATATYPE>
91
92 template<typename DATATYPE>
93 int AMPI_Reduce_global(MEDI_OPTIONAL_CONST typename DATATYPE::Type* sendbuf, typename DATATYPE::Type* recvbuf, int count, DATATYPE* datatype, AMPI_Op op, int root, AMPI_Comm comm);
94 template<typename SENDTYPE, typename RECVTYPE>
95 int AMPI_Gather(MEDI_OPTIONAL_CONST typename SENDTYPE::Type* sendbuf, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::Type* recvbuf, int recvcount, RECVTYPE* recvtype, int root, AMPI_Comm comm);
96 template<typename DATATYPE>
97 int AMPI_Allreduce_global(MEDI_OPTIONAL_CONST typename DATATYPE::Type* sendbuf, typename DATATYPE::Type* recvbuf, int count, DATATYPE* datatype, AMPI_Op op, AMPI_Comm comm);
98 template<typename SENDTYPE, typename RECVTYPE>
99 int AMPI_Allgather(MEDI_OPTIONAL_CONST typename SENDTYPE::Type* sendbuf, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::Type* recvbuf, int recvcount, RECVTYPE* recvtype, AMPI_Comm comm);
100 template<typename DATATYPE>
101 int AMPI_Ireduce_global(MEDI_OPTIONAL_CONST typename DATATYPE::Type* sendbuf, typename DATATYPE::Type* recvbuf, int count, DATATYPE* datatype, AMPI_Op op, int root, AMPI_Comm comm, AMPI_Request* request);
102 template<typename SENDTYPE, typename RECVTYPE>
103 int AMPI_Igather(MEDI_OPTIONAL_CONST typename SENDTYPE::Type* sendbuf, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::Type* recvbuf, int recvcount, RECVTYPE* recvtype, int root, AMPI_Comm comm, AMPI_Request* request);
104 template<typename DATATYPE>
105 int AMPI_Iallreduce_global(MEDI_OPTIONAL_CONST typename DATATYPE::Type* sendbuf, typename DATATYPE::Type* recvbuf, int count, DATATYPE* datatype, AMPI_Op op, AMPI_Comm comm, AMPI_Request* request);
106 template<typename SENDTYPE, typename RECVTYPE>
107 int AMPI_Iallgather(MEDI_OPTIONAL_CONST typename SENDTYPE::Type* sendbuf, int sendcount, SENDTYPE* sendtype, typename RECVTYPE::Type* recvbuf, int recvcount, RECVTYPE* recvtype, AMPI_Comm comm, AMPI_Request* request);
108
109 template<typename DATATYPE>
110 inline void performReduce(typename DATATYPE::Type* tempbuf, typename DATATYPE::Type* recvbuf, int count, DATATYPE* datatype, AMPI_Op op, int root, AMPI_Comm comm, int reduceSize) {
111 int commRank = getCommRank(comm);
112
113 if(-1 == root || root == commRank) {
114 datatype->performReduce(tempbuf, recvbuf, count, op, reduceSize);
115 }
116 }
117
118 template<typename DATATYPE>
119 inline int GatherAndPerformOperationLocal(MEDI_OPTIONAL_CONST typename DATATYPE::Type* sendbuf, typename DATATYPE::Type* recvbuf, int count, DATATYPE* datatype, AMPI_Op op, int root, AMPI_Comm comm, int reduceSize) {
120 int commSize = getCommSize(comm);
121 int commRank = getCommRank(comm);
122
123 typename DATATYPE::Type* tempbuf = NULL;
124 if(-1 == root || root == commRank) {
125 datatype->createTypeBuffer(tempbuf, count * commSize);
126 }
127
128 MEDI_OPTIONAL_CONST typename DATATYPE::Type* sendbufGather = sendbuf;
129 if(AMPI_IN_PLACE == sendbuf) {
130 sendbufGather = recvbuf;
131 }
132
133 int rValue;
134 if(-1 == root) {
135 rValue = AMPI_Allgather<DATATYPE, DATATYPE>(sendbufGather, count, datatype, tempbuf, count, datatype, comm);
136 } else {
137 rValue = AMPI_Gather<DATATYPE, DATATYPE>(sendbufGather, count, datatype, tempbuf, count, datatype, root, comm);
138 }
139
140 performReduce<DATATYPE>(tempbuf, recvbuf, count, datatype, op, root, comm, reduceSize);
141
142 if(-1 == root || root == commRank) {
143 datatype->deleteTypeBuffer(tempbuf, count * commSize);
144 }
145
146 return rValue;
147 }
148
149 template<typename DATATYPE>
151
153 // first call the orignal function
154 h->origFunc(h->origHandle);
155
156 // call custom function if present
157 if (h->customOpPriorToLocalReduce != nullptr) {
159 }
160
161 int commSize = getCommSize(h->comm);
162 int commRank = getCommRank(h->comm);
163
164 performReduce<DATATYPE>(h->tempbuf, h->recvbuf, h->count, h->datatype, h->op, h->root, h->comm, h->reduceSize);
165
166 if(-1 == h->root || h->root == commRank) {
167 h->datatype->deleteTypeBuffer(h->tempbuf, h->count * commSize);
168 }
169
170 delete h;
171
172 return 0;
173 }
174
175 template<typename DATATYPE>
176 inline int IgatherAndPerformOperationLocal(MEDI_OPTIONAL_CONST typename DATATYPE::Type* sendbuf, typename DATATYPE::Type* recvbuf, int count, DATATYPE* datatype, AMPI_Op op, int root, AMPI_Comm comm, AMPI_Request* request, int reduceSize) {
177 int commSize = getCommSize(comm);
178 int commRank = getCommRank(comm);
179
180 typename DATATYPE::Type* tempbuf = NULL;
181 if(-1 == root || root == commRank) {
182 datatype->createTypeBuffer(tempbuf, count * commSize);
183 }
184
185 MEDI_OPTIONAL_CONST typename DATATYPE::Type* sendbufGather = sendbuf;
186 if(AMPI_IN_PLACE == sendbuf) {
187 sendbufGather = recvbuf;
188 }
189
190 int rValue;
191 if(-1 == root) {
192 rValue = AMPI_Iallgather<DATATYPE, DATATYPE>(sendbufGather, count, datatype, tempbuf, count, datatype, comm, request);
193 } else {
194 rValue = AMPI_Igather<DATATYPE, DATATYPE>(sendbufGather, count, datatype, tempbuf, count, datatype, root, comm, request);
195 }
197 curHandle->comm = comm;
198 curHandle->root = root;
199 curHandle->count = count;
200 curHandle->op = op;
201 curHandle->datatype = datatype;
202 curHandle->recvbuf = recvbuf;
203 curHandle->tempbuf = tempbuf;
204 curHandle->origHandle = request->handle;
205 curHandle->origFunc = request->func;
206 curHandle->reduceSize = reduceSize;
207 curHandle->toolHandle = request->handle->toolHandle;
208
209 // set our own handle now to the request
210 request->handle = curHandle;
211 request->func = (ContinueFunction)IgatherAndPerformOperationLocal_finish<DATATYPE>;
212
213 return rValue;
214 }
215
216 template<typename DATATYPE>
217 inline bool addCustomOperationPriorToLocalReduce(AMPI_Request* request, CustomFunction func, void* data) {
218
220
221 if (reduceHandle != nullptr) {
222 reduceHandle->customOpPriorToLocalReduce = func;
223 reduceHandle->customDataPriorToLocalReduce = data;
224 return true;
225 } else { // no local reduce performed
226 return false;
227 }
228 }
229
230 template<typename DATATYPE>
241
242 template<typename DATATYPE>
244
246 // first call the orignal function
247 h->origFunc(h->origHandle);
248
249 if(h->root != getCommRank(h->comm)) {
250 h->datatype->deleteTypeBuffer(h->tempBuf, h->count);
251 }
252
253 delete h;
254
255 return 0;
256 }
257
258 template<typename DATATYPE>
259 inline int AMPI_Reduce(MEDI_OPTIONAL_CONST typename DATATYPE::Type* sendbuf, typename DATATYPE::Type* recvbuf, int count, DATATYPE* datatype, AMPI_Op op, int root, AMPI_Comm comm) {
260 AMPI_Op convOp = datatype->getADTool().convertOperator(op);
261
262 if(!datatype->getADTool().isActiveType()) {
263 return AMPI_Reduce_global<DATATYPE>(sendbuf, recvbuf, count, datatype, convOp, root, comm);
264 } else if(convOp.hasAdjoint) {
265 // operator has an adjoint formulation
266 if(convOp.requiresPrimalSend) {
267 // need to modify the call to an allreduce
268 typename DATATYPE::Type* tempBuf = recvbuf;
269 if(root != getCommRank(comm)) {
270 datatype->createTypeBuffer(tempBuf, count);
271 }
272 int result = AMPI_Allreduce_global<DATATYPE>(sendbuf, tempBuf, count, datatype, convOp, comm);
273 if(root != getCommRank(comm)) {
274 datatype->deleteTypeBuffer(tempBuf, count);
275 }
276
277 return result;
278 } else {
279 // just perfrom the normal call
280 return AMPI_Reduce_global<DATATYPE>(sendbuf, recvbuf, count, datatype, convOp, root, comm);
281 }
282 } else {
283 // perform a gather and apply the operator locally
284 return GatherAndPerformOperationLocal(sendbuf, recvbuf, count, datatype, convOp, root, comm, getCommSize(comm));
285 }
286 }
287
288 template<typename DATATYPE>
289 inline int AMPI_Ireduce(MEDI_OPTIONAL_CONST typename DATATYPE::Type* sendbuf, typename DATATYPE::Type* recvbuf, int count, DATATYPE* datatype, AMPI_Op op, int root, AMPI_Comm comm, AMPI_Request* request) {
290 AMPI_Op convOp = datatype->getADTool().convertOperator(op);
291
292 if(!datatype->getADTool().isActiveType()) {
293 return AMPI_Ireduce_global<DATATYPE>(sendbuf, recvbuf, count, datatype, convOp, root, comm, request);
294 } else if(convOp.hasAdjoint) {
295 if(convOp.requiresPrimalSend) {
296 // need to modify the call to an allreduce
297 typename DATATYPE::Type* tempBuf = recvbuf;
298 if(root != getCommRank(comm)) {
299 datatype->createTypeBuffer(tempBuf, count);
300 }
301 int result = AMPI_Iallreduce_global<DATATYPE>(sendbuf, tempBuf, count, datatype, convOp, comm, request);
303 curHandle->tempBuf = tempBuf;
304 curHandle->comm = comm;
305 curHandle->root = root;
306 curHandle->count = count;
307 curHandle->datatype = datatype;
308 curHandle->origHandle = request->handle;
309 curHandle->origFunc = request->func;
310 curHandle->toolHandle = request->handle->toolHandle;
311
312 // set our own handle now to the request
313 request->handle = curHandle;
314 request->func = (ContinueFunction)AMPI_Ireduce_modified_finish<DATATYPE>;
315
316 return result;
317 } else {
318 // just perfrom the normal call
319 return AMPI_Ireduce_global<DATATYPE>(sendbuf, recvbuf, count, datatype, convOp, root, comm, request);
320 }
321 } else {
322 // perform a gather and apply the operator locally
323 return IgatherAndPerformOperationLocal(sendbuf, recvbuf, count, datatype, convOp, root, comm, request, getCommSize(comm));
324 }
325 }
326
327 template<typename DATATYPE>
328 inline int AMPI_Allreduce(MEDI_OPTIONAL_CONST typename DATATYPE::Type* sendbuf, typename DATATYPE::Type* recvbuf, int count, DATATYPE* datatype, AMPI_Op op, AMPI_Comm comm) {
329 AMPI_Op convOp = datatype->getADTool().convertOperator(op);
330
331 if(convOp.hasAdjoint || !datatype->getADTool().isActiveType()) {
332 return AMPI_Allreduce_global<DATATYPE>(sendbuf, recvbuf, count, datatype, convOp, comm);
333 } else {
334 // perform a gather and apply the operator locally
335 return GatherAndPerformOperationLocal(sendbuf, recvbuf, count, datatype, convOp, -1, comm, getCommSize(comm));
336 }
337 }
338
339 template<typename DATATYPE>
340 inline int AMPI_Iallreduce(MEDI_OPTIONAL_CONST typename DATATYPE::Type* sendbuf, typename DATATYPE::Type* recvbuf, int count, DATATYPE* datatype, AMPI_Op op, AMPI_Comm comm, AMPI_Request* request) {
341 AMPI_Op convOp = datatype->getADTool().convertOperator(op);
342
343 if(convOp.hasAdjoint || !datatype->getADTool().isActiveType()) {
344 return AMPI_Iallreduce_global<DATATYPE>(sendbuf, recvbuf, count, datatype, convOp, comm, request);
345 } else {
346 // perform a gather and apply the operator locally
347 return IgatherAndPerformOperationLocal(sendbuf, recvbuf, count, datatype, convOp, -1, comm, request, getCommSize(comm));
348 }
349 }
350
351 template<typename DATATYPE>
352 inline int AMPI_Exscan(const typename DATATYPE::Type* sendbuf, typename DATATYPE::Type* recvbuf, int count, DATATYPE* datatype, AMPI_Op op, AMPI_Comm comm) {
353 AMPI_Op convOp = datatype->getADTool().convertOperator(op);
354
355 if(!datatype->getADTool().isActiveType()) {
356 return MPI_Exscan(sendbuf, recvbuf, count, datatype->getMpiType(), convOp.primalFunction, comm);
357 } else {
358 // perform a gather and apply the operator locally
359 return GatherAndPerformOperationLocal(sendbuf, recvbuf, count, datatype, convOp, -1, comm, getCommRank(comm));
360 }
361 }
362
363 template<typename DATATYPE>
364 inline int AMPI_Iexscan(const typename DATATYPE::Type* sendbuf, typename DATATYPE::Type* recvbuf, int count, DATATYPE* datatype, AMPI_Op op, AMPI_Comm comm, AMPI_Request* request) {
365 AMPI_Op convOp = datatype->getADTool().convertOperator(op);
366
367 if(!datatype->getADTool().isActiveType()) {
368 return MPI_Iexscan(sendbuf, recvbuf, count, datatype->getMpiType(), convOp.primalFunction, comm, &request->request);
369 } else {
370 // perform a gather and apply the operator locally
371 return IgatherAndPerformOperationLocal(sendbuf, recvbuf, count, datatype, convOp, -1, comm, request, getCommRank(comm));
372 }
373 }
374
375 template<typename DATATYPE>
376 inline int AMPI_Scan(const typename DATATYPE::Type* sendbuf, typename DATATYPE::Type* recvbuf, int count, DATATYPE* datatype, AMPI_Op op, AMPI_Comm comm) {
377 AMPI_Op convOp = datatype->getADTool().convertOperator(op);
378
379 if(!datatype->getADTool().isActiveType()) {
380 return MPI_Scan(sendbuf, recvbuf, count, datatype->getMpiType(), convOp.primalFunction, comm);
381 } else {
382 // perform a gather and apply the operator locally
383 return GatherAndPerformOperationLocal(sendbuf, recvbuf, count, datatype, convOp, -1, comm, getCommRank(comm) + 1);
384 }
385 }
386
387 template<typename DATATYPE>
388 inline int AMPI_Iscan(const typename DATATYPE::Type* sendbuf, typename DATATYPE::Type* recvbuf, int count, DATATYPE* datatype, AMPI_Op op, AMPI_Comm comm, AMPI_Request* request) {
389 AMPI_Op convOp = datatype->getADTool().convertOperator(op);
390
391 if(!datatype->getADTool().isActiveType()) {
392 return MPI_Iscan(sendbuf, recvbuf, count, datatype->getMpiType(), convOp.primalFunction, comm, &request->request);
393 } else {
394 // perform a gather and apply the operator locally
395 return IgatherAndPerformOperationLocal(sendbuf, recvbuf, count, datatype, convOp, -1, comm, request, getCommRank(comm) + 1);
396 }
397 }
398
399 inline void AMPI_Init_common() {
400 initTypes();
402 }
403
404
405 inline int AMPI_Init(int* argc, char*** argv) {
406 int result = MPI_Init(argc, argv);
407
409
410 return result;
411 }
412
413 inline int AMPI_Init_thread(int *argc, char ***argv, int required, int *provided) {
414 int result = MPI_Init_thread(argc, argv, required, provided);
415
417
418 return result;
419 }
420
421 template<typename DATATYPE>
422 inline int AMPI_Reduce_local(MEDI_OPTIONAL_CONST typename DATATYPE::Type* inbuf, typename DATATYPE::Type* inoutbuf, int count, DATATYPE* datatype, AMPI_Comm comm, AMPI_Op op) {
423 AMPI_Op convOp = datatype->getADTool().convertOperator(op);
424
425 return MPI_Reduce_local(inbuf, inoutbuf, count, datatype->getMpiType(), convOp.primalFunction);
426 }
427
428
429}
#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_Init_common()
Definition wrappers.hpp:399
int AMPI_Reduce_global(const typename DATATYPE::Type *sendbuf, typename DATATYPE::Type *recvbuf, int count, DATATYPE *datatype, AMPI_Op op, int root, MPI_Comm comm)
Definition ampiFunctions.hpp:10845
int AMPI_Allgather(const typename SENDTYPE::Type *sendbuf, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::Type *recvbuf, int recvcount, RECVTYPE *recvtype, MPI_Comm comm)
Definition ampiFunctions.hpp:4061
void performReduce(typename DATATYPE::Type *tempbuf, typename DATATYPE::Type *recvbuf, int count, DATATYPE *datatype, AMPI_Op op, int root, MPI_Comm comm, int reduceSize)
Definition wrappers.hpp:110
int MPI_Ibcast_wrap(void *bufferSend, void *bufferRecv, int count, MPI_Datatype type, int root, MPI_Comm comm, MPI_Request *request)
Definition wrappers.hpp:68
int AMPI_Iexscan(const typename DATATYPE::Type *sendbuf, typename DATATYPE::Type *recvbuf, int count, DATATYPE *datatype, AMPI_Op op, MPI_Comm comm, AMPI_Request *request)
Definition wrappers.hpp:364
void(* CustomFunction)(void *data)
Definition typeDefinitions.h:51
int AMPI_Gather(const typename SENDTYPE::Type *sendbuf, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::Type *recvbuf, int recvcount, RECVTYPE *recvtype, int root, MPI_Comm comm)
Definition ampiFunctions.hpp:5758
void initializeOperators()
Definition ampiDefinitions.cpp:196
int AMPI_Reduce_local(const typename DATATYPE::Type *inbuf, typename DATATYPE::Type *inoutbuf, int count, DATATYPE *datatype, MPI_Comm comm, AMPI_Op op)
Definition wrappers.hpp:422
int AMPI_Igather(const typename SENDTYPE::Type *sendbuf, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::Type *recvbuf, int recvcount, RECVTYPE *recvtype, int root, MPI_Comm comm, AMPI_Request *request)
Definition ampiFunctions.hpp:8753
int AMPI_Ibcast_wrap(typename DATATYPE::Type *bufferSend, typename DATATYPE::Type *bufferRecv, int count, DATATYPE *datatype, int root, MPI_Comm comm, AMPI_Request *request)
Definition ampiFunctions.hpp:8375
int AMPI_Exscan(const typename DATATYPE::Type *sendbuf, typename DATATYPE::Type *recvbuf, int count, DATATYPE *datatype, AMPI_Op op, MPI_Comm comm)
Definition wrappers.hpp:352
const AMPI_IN_PLACE_IMPL AMPI_IN_PLACE
This structure is able to be cast to any type.
Definition ampi.cpp:40
int MPI_Bcast_wrap(void *bufferSend, void *bufferRecv, int count, MPI_Datatype type, int root, MPI_Comm comm)
Definition wrappers.hpp:52
int AMPI_Iallgather(const typename SENDTYPE::Type *sendbuf, int sendcount, SENDTYPE *sendtype, typename RECVTYPE::Type *recvbuf, int recvcount, RECVTYPE *recvtype, MPI_Comm comm, AMPI_Request *request)
Definition ampiFunctions.hpp:6410
int IgatherAndPerformOperationLocal(const typename DATATYPE::Type *sendbuf, typename DATATYPE::Type *recvbuf, int count, DATATYPE *datatype, AMPI_Op op, int root, MPI_Comm comm, AMPI_Request *request, int reduceSize)
Definition wrappers.hpp:176
int GatherAndPerformOperationLocal(const typename DATATYPE::Type *sendbuf, typename DATATYPE::Type *recvbuf, int count, DATATYPE *datatype, AMPI_Op op, int root, MPI_Comm comm, int reduceSize)
Definition wrappers.hpp:119
int AMPI_Init_thread(int *argc, char ***argv, int required, int *provided)
Definition wrappers.hpp:413
int AMPI_Allreduce(const typename DATATYPE::Type *sendbuf, typename DATATYPE::Type *recvbuf, int count, DATATYPE *datatype, AMPI_Op op, MPI_Comm comm)
Definition wrappers.hpp:328
int AMPI_Ireduce(const typename DATATYPE::Type *sendbuf, typename DATATYPE::Type *recvbuf, int count, DATATYPE *datatype, AMPI_Op op, int root, MPI_Comm comm, AMPI_Request *request)
Definition wrappers.hpp:289
int IgatherAndPerformOperationLocal_finish(HandleBase *handle)
Definition wrappers.hpp:150
int AMPI_Allreduce_global(const typename DATATYPE::Type *sendbuf, typename DATATYPE::Type *recvbuf, int count, DATATYPE *datatype, AMPI_Op op, MPI_Comm comm)
Definition ampiFunctions.hpp:4623
void initTypes()
Definition ampiDefinitions.cpp:241
int getCommSize(MPI_Comm comm)
Helper function that gets the number of ranks from the communicator.
Definition mpiTools.h:64
int getCommRank(MPI_Comm comm)
Helper function that gets the own rank number from the communicator.
Definition mpiTools.h:52
int AMPI_Reduce(const typename DATATYPE::Type *sendbuf, typename DATATYPE::Type *recvbuf, int count, DATATYPE *datatype, AMPI_Op op, int root, MPI_Comm comm)
Definition wrappers.hpp:259
int AMPI_Ibcast(typename DATATYPE::Type *buffer, int count, DATATYPE *datatype, int root, MPI_Comm comm, AMPI_Request *request)
Definition wrappers.hpp:64
int AMPI_Bcast_wrap(typename DATATYPE::Type *bufferSend, typename DATATYPE::Type *bufferRecv, int count, DATATYPE *datatype, int root, MPI_Comm comm)
Definition ampiFunctions.hpp:5483
bool addCustomOperationPriorToLocalReduce(AMPI_Request *request, CustomFunction func, void *data)
Definition wrappers.hpp:217
int AMPI_Iallreduce_global(const typename DATATYPE::Type *sendbuf, typename DATATYPE::Type *recvbuf, int count, DATATYPE *datatype, AMPI_Op op, MPI_Comm comm, AMPI_Request *request)
Definition ampiFunctions.hpp:7196
int AMPI_Iscan(const typename DATATYPE::Type *sendbuf, typename DATATYPE::Type *recvbuf, int count, DATATYPE *datatype, AMPI_Op op, MPI_Comm comm, AMPI_Request *request)
Definition wrappers.hpp:388
int(* ContinueFunction)(HandleBase *h)
Definition typeDefinitions.h:48
int AMPI_Bcast(typename DATATYPE::Type *buffer, int count, DATATYPE *datatype, int root, MPI_Comm comm)
Definition wrappers.hpp:48
int AMPI_Init(int *argc, char ***argv)
Definition wrappers.hpp:405
int AMPI_Ireduce_global(const typename DATATYPE::Type *sendbuf, typename DATATYPE::Type *recvbuf, int count, DATATYPE *datatype, AMPI_Op op, int root, MPI_Comm comm, AMPI_Request *request)
Definition ampiFunctions.hpp:9611
int AMPI_Scan(const typename DATATYPE::Type *sendbuf, typename DATATYPE::Type *recvbuf, int count, DATATYPE *datatype, AMPI_Op op, MPI_Comm comm)
Definition wrappers.hpp:376
int AMPI_Iallreduce(const typename DATATYPE::Type *sendbuf, typename DATATYPE::Type *recvbuf, int count, DATATYPE *datatype, AMPI_Op op, MPI_Comm comm, AMPI_Request *request)
Definition wrappers.hpp:340
int AMPI_Ireduce_modified_finish(HandleBase *handle)
Definition wrappers.hpp:243
Definition wrappers.hpp:75
DATATYPE::Type * tempbuf
Definition wrappers.hpp:82
ContinueFunction origFunc
Definition wrappers.hpp:85
CustomFunction customOpPriorToLocalReduce
Definition wrappers.hpp:88
HandleBase * origHandle
Definition wrappers.hpp:84
DATATYPE * datatype
Definition wrappers.hpp:79
int reduceSize
Definition wrappers.hpp:86
DATATYPE::Type * recvbuf
Definition wrappers.hpp:81
AMPI_Op op
Definition wrappers.hpp:80
void * customDataPriorToLocalReduce
Definition wrappers.hpp:89
MPI_Comm comm
Definition wrappers.hpp:76
int root
Definition wrappers.hpp:77
int count
Definition wrappers.hpp:78
Definition wrappers.hpp:231
ContinueFunction origFunc
Definition wrappers.hpp:239
DATATYPE * datatype
Definition wrappers.hpp:233
DATATYPE::Type * tempBuf
Definition wrappers.hpp:232
int root
Definition wrappers.hpp:234
MPI_Comm comm
Definition wrappers.hpp:236
int count
Definition wrappers.hpp:235
HandleBase * origHandle
Definition wrappers.hpp:238
Structure for the special handling of the MPI_Op structure.
Definition op.hpp:50
MPI_Op primalFunction
The mpi operator for the unmodified AD types. The AD tool needs to record all operations that are eva...
Definition op.hpp:69
bool requiresPrimalSend
Indicates if the primal on the receiving side needs to be send to the sending side.
Definition op.hpp:60
bool hasAdjoint
Indicates if the user has provided a specialized adjoint handling for the operator.
Definition op.hpp:92
Definition async.hpp:45
AsyncHandle * handle
Definition async.hpp:47
MPI_Request request
Definition async.hpp:46
ContinueFunction func
Definition async.hpp:48
Definition async.hpp:101
AsyncAdjointHandle * toolHandle
Definition async.hpp:103
Definition typeDefinitions.h:53