MeDiPack  1.3.0
A Message Differentiation Package
SciComp TU Kaiserslautern
Loading...
Searching...
No Matches
constructedDatatypes.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 <cstdlib>
32
33#include "../macros.h"
34#include "typeInterface.hpp"
35#include "../exceptions.hpp"
36
40namespace medi {
41
48 class MpiStructType final : public MpiTypeInterface {
49
50 private:
51 bool modificationRequired;
52 int valuesPerElement;
53
54 const ADToolInterface* adInterface;
55
56 size_t typeExtent;
57 size_t typeOffset;
58 size_t modifiedExtent;
59
60 int nTypes;
61 int* blockLengths;
62 int* blockOffsets;
63 int* modifiedBlockOffsets;
64 MpiTypeInterface** types;
65
66
67 public:
68 typedef void Type;
69 typedef void ModifiedType;
70 typedef void AdjointType;
71 typedef void PrimalType;
72 typedef void IndexType;
73
74 private:
75 void cloneInternal(const MpiStructType* other) {
76 nTypes = other->nTypes;
77 blockLengths = new int [nTypes];
78 blockOffsets = new int [nTypes];
79 types = new MpiTypeInterface* [nTypes];
80
81 if(nullptr != other->modifiedBlockOffsets) {
82 modifiedBlockOffsets = new int [nTypes];
83 }
84
85 for(int i = 0; i < nTypes; ++i) {
86 blockLengths[i] = other->blockLengths[i];
87 blockOffsets[i] = other->blockOffsets[i];
88 types[i] = other->types[i]->clone();
89
90 if(nullptr != modifiedBlockOffsets) {
91 modifiedBlockOffsets[i] = other->modifiedBlockOffsets[i];
92 }
93 }
94 }
95
96 public:
97
98
100 MpiTypeInterface(MPI_INT, MPI_INT)
101 {
102 adInterface = other->adInterface;
103 typeExtent = other->typeExtent;
104 typeOffset = other->typeOffset;
105 modifiedExtent = other->modifiedExtent;
106
107 cloneInternal(other);
108
109 MPI_Datatype type;
110 MPI_Datatype modType;
111
112 MPI_Type_dup(this->getMpiType(), &type);
113 if(this->getMpiType() != this->getModifiedMpiType()) {
114 MPI_Type_dup(this->getModifiedMpiType(), &modType);
115 } else {
116 modType = type;
117 }
118
119 setMpiTypes(type, modType);
120
121 }
122
123 MpiStructType(const MpiStructType* other, size_t offset, size_t extent) :
124 MpiTypeInterface(MPI_INT, MPI_INT)
125 {
126 adInterface = other->adInterface;
127 typeExtent = extent;
128 typeOffset = offset;
129 if(nullptr != modifiedBlockOffsets) {
130 modifiedExtent = other->modifiedExtent;
131 } else {
132 modifiedExtent = extent;
133 }
134
135 cloneInternal(other);
136
137 MPI_Datatype type;
138 MPI_Datatype modType;
139
140 MPI_Type_create_resized(this->getMpiType(), offset, extent, &type);
141 if(this->getMpiType() != this->getModifiedMpiType()) {
142 MPI_Type_dup(this->getModifiedMpiType(), &modType);
143 } else {
144 modType = type;
145 }
146
147 setMpiTypes(type, modType);
148 }
149
150
151 MpiStructType(int count, MEDI_OPTIONAL_CONST int* array_of_blocklengths, MEDI_OPTIONAL_CONST MPI_Aint* array_of_displacements, MpiTypeInterface* const * array_of_types) :
152 MpiTypeInterface(MPI_INT, MPI_INT) {
153
154 MPI_Datatype* mpiTypes = new MPI_Datatype[count];
155
156 MPI_Datatype newMpiType;
157 MPI_Datatype newModMpiType;
158
159 nTypes = count;
160 blockLengths = new int[count];
161 blockOffsets = new int[count];
162 types = new MpiTypeInterface*[count];
163
164 // check if a modified buffer is required and populate the mpiTypes as well as the arrayes
165 modificationRequired = false;
166 valuesPerElement = 0;
167 adInterface = &array_of_types[0]->getADTool(); // Set AD interface of first type as the default.
168 for(int i = 0; i < count; ++i) {
169 modificationRequired |= array_of_types[i]->isModifiedBufferRequired();
170 if(array_of_types[i]->getADTool().isActiveType()) {
171 adInterface = &array_of_types[i]->getADTool();
172 }
173
174 valuesPerElement += array_of_types[i]->computeActiveElements(array_of_blocklengths[i]);
175
176 blockLengths[i] = array_of_blocklengths[i];
177 blockOffsets[i] = array_of_displacements[i];
178 mpiTypes[i] = array_of_types[i]->getMpiType();
179 types[i] = array_of_types[i]->clone();
180 }
181
182 MPI_Type_create_struct(count, array_of_blocklengths, array_of_displacements, mpiTypes, &newMpiType);
183
184 if(modificationRequired) {
185 MPI_Aint* modifiedDisplacements = new MPI_Aint[count + 1]; // We might need to add padding so add an extra element
186 MPI_Datatype* modifiedMpiTypes = new MPI_Datatype[count + 1]; // We might need to add padding so add an extra element
187 int* modifiedArrayLength = new int[count + 1]; // We might need to add padding so add an extra element
188
189 modifiedBlockOffsets = new int[count];
190
191 MPI_Aint totalDisplacement = 0;
192 for(int i = 0; i < count; ++i) {
193 MPI_Aint curLowerBound;
194 MPI_Aint curExtent;
195
196 modifiedMpiTypes[i] = array_of_types[i]->getModifiedMpiType();
197
198#if MEDI_MPI_TARGET < MEDI_MPI_VERSION_2_0
199 MPI_Type_lb(modifiedMpiTypes[i], &curLowerBound);
200 MPI_Type_extent(modifiedMpiTypes[i], &curExtent);
201#else
202 MPI_Type_get_extent(modifiedMpiTypes[i], &curLowerBound, &curExtent);
203#endif
204
205 mediAssert(0 == curLowerBound); // The modified types are always packed without any holes.
206 modifiedDisplacements[i] = totalDisplacement;
207 modifiedBlockOffsets[i] = totalDisplacement;
208 modifiedArrayLength[i] = array_of_blocklengths[i];
209 totalDisplacement += curExtent * array_of_blocklengths[i];
210 }
211
212 // TODO: This code assumes a 64-bit machine and if the last member of a struct is a byte, that these are
213 // padding bytes
214 size_t paddingBytes = totalDisplacement % sizeof(double);
215 int numberOfTypes = count; // This is obvious default we might increase this in the ifs below
216 if(paddingBytes != 0) {
217 // yuck we have to add padding
218
219 // check if the last member was a byte, if yes recompute the total displacement
220 if(modifiedMpiTypes[count - 1] == MPI_BYTE) {
221 totalDisplacement -= array_of_blocklengths[count - 1]; // extent is one
222 paddingBytes = totalDisplacement % sizeof(double);
223
224 if(paddingBytes != 0) {
225 // we still have padding bytes so change the array size
226 modifiedArrayLength[count - 1] = paddingBytes;
227 totalDisplacement += paddingBytes;
228 } else {
229 // no padding bytes left remove the padding bytes from the creation
230 modifiedArrayLength[count - 1] = 0;
231 }
232
233 // since we assume now that the last type are padding bytes, we can remove it from the iterations
234 blockLengths[count - 1] = 0;
235 } else {
236 // The last type is now padding type so we have to add a new type to the lists
237 modifiedMpiTypes[count] = MPI_BYTE;
238 modifiedDisplacements[count] = totalDisplacement;
239 modifiedArrayLength[count] = paddingBytes;
240 totalDisplacement += paddingBytes;
241
242 numberOfTypes += 1;
243 }
244 }
245
246 MPI_Type_create_struct(numberOfTypes, modifiedArrayLength, modifiedDisplacements, modifiedMpiTypes, &newModMpiType);
247
248 delete [] modifiedMpiTypes;
249 delete [] modifiedDisplacements;
250 delete [] modifiedArrayLength;
251 } else {
252
253 modifiedBlockOffsets = new int[count];
254 for(int i = 0; i < count; ++i) {
255 modifiedBlockOffsets[i] = blockOffsets[i];
256 }
257 newModMpiType = newMpiType;
258 }
259
260 MPI_Aint lb = 0;
261 MPI_Aint ext = 0;
262#if MEDI_MPI_TARGET < MEDI_MPI_VERSION_2_0
263 MPI_Type_lb(newMpiType, &lb);
264 MPI_Type_extent(newMpiType, &ext);
265#else
266 MPI_Type_get_extent(newMpiType, &lb, &ext);
267#endif
268 typeOffset = lb;
269 typeExtent = ext;
270
271#if MEDI_MPI_TARGET < MEDI_MPI_VERSION_2_0
272 MPI_Type_lb(newModMpiType, &lb);
273 MPI_Type_extent(newModMpiType, &ext);
274#else
275 MPI_Type_get_extent(newModMpiType, &lb, &ext);
276#endif
277 modifiedExtent = ext;
278
279 setMpiTypes(newMpiType, newModMpiType);
280
281 delete [] mpiTypes;
282 }
283
285
286 MPI_Datatype temp;
287 if(this->getModifiedMpiType() != this->getMpiType()) {
288 temp = this->getModifiedMpiType();
289 MPI_Type_free(&temp);
290 }
291 temp = this->getMpiType();
292 MPI_Type_free(&temp);
293
294 for(int i = 0; i < nTypes; ++i) {
295 delete types[i];
296 }
297
298 if(nullptr != modifiedBlockOffsets) {
299 delete [] modifiedBlockOffsets;
300 }
301 delete [] blockOffsets;
302 delete [] blockLengths;
303 delete [] types;
304 }
305
306 int computeBufOffset(size_t element) const {
307 return (int)(element * typeExtent);
308 }
309
310 int computeModOffset(size_t element) const {
311 return (int)(element * modifiedExtent);
312 }
313
314 const void* computeBufferPointer(const void* buf, size_t offset) const {
315 return (const void*)((const char*) buf + offset);
316 }
317
318 void* computeBufferPointer(void* buf, size_t offset) const {
319 return (void*)((char*) buf + offset);
320 }
321
323 return modificationRequired;
324 }
325
326 int computeActiveElements(const int count) const {
327 return count * valuesPerElement;
328 }
329
330
331 const ADToolInterface& getADTool() const {
332 return *adInterface;
333 }
334
335 void copyIntoModifiedBuffer(const void* buf, size_t bufOffset, void* bufMod, size_t bufModOffset, int elements) const {
336 for(int i = 0; i < elements; ++i) {
337 int totalBufOffset = computeBufOffset(i + bufOffset);
338 int totalModOffset = computeModOffset(i + bufModOffset);
339
340 for(int curType = 0; curType < nTypes; ++curType) {
341 types[curType]->copyIntoModifiedBuffer(computeBufferPointer(buf, totalBufOffset + blockOffsets[curType]), 0, computeBufferPointer(bufMod, totalModOffset + modifiedBlockOffsets[curType]), 0, blockLengths[curType]);
342 }
343 }
344 }
345
346 void copyFromModifiedBuffer(void* buf, size_t bufOffset, const void* bufMod, size_t bufModOffset, int elements) const {
347 for(int i = 0; i < elements; ++i) {
348 int totalBufOffset = computeBufOffset(i + bufOffset);
349 int totalModOffset = computeModOffset(i + bufModOffset);
350
351
352 for(int curType = 0; curType < nTypes; ++curType) {
353 types[curType]->copyFromModifiedBuffer(computeBufferPointer(buf, totalBufOffset + blockOffsets[curType]), 0, computeBufferPointer(bufMod, totalModOffset + modifiedBlockOffsets[curType]), 0, blockLengths[curType]);
354 }
355 }
356 }
357
358 void getIndices(const void* buf, size_t bufOffset, void* indices, size_t bufModOffset, int elements) const {
359 int totalIndexOffset = computeActiveElements(bufModOffset); // indices are lineralized and counted up in the loop
360
361 for(int i = 0; i < elements; ++i) {
362 int totalBufOffset = computeBufOffset(i + bufOffset);
363
364 for(int curType = 0; curType < nTypes; ++curType) {
365 types[curType]->getIndices(computeBufferPointer(buf, totalBufOffset + blockOffsets[curType]), 0, indices, totalIndexOffset, blockLengths[curType]);
366 totalIndexOffset += types[curType]->computeActiveElements(blockLengths[curType]);
367 }
368 }
369 }
370
371 void registerValue(void* buf, size_t bufOffset, void* indices, void* oldPrimals, size_t bufModOffset, int elements) const {
372 int totalIndexOffset = computeActiveElements(bufModOffset); // indices are lineralized and counted up in the loop
373
374 for(int i = 0; i < elements; ++i) {
375 int totalBufOffset = computeBufOffset(i + bufOffset);
376
377 for(int curType = 0; curType < nTypes; ++curType) {
378 types[curType]->registerValue(computeBufferPointer(buf, totalBufOffset + blockOffsets[curType]), 0, indices, oldPrimals, totalIndexOffset, blockLengths[curType]);
379 totalIndexOffset += types[curType]->computeActiveElements(blockLengths[curType]);
380 }
381 }
382 }
383
384 void clearIndices(void* buf, size_t bufOffset, int elements) const {
385 for(int i = 0; i < elements; ++i) {
386 int totalBufOffset = computeBufOffset(i + bufOffset);
387
388 for(int curType = 0; curType < nTypes; ++curType) {
389 types[curType]->clearIndices(computeBufferPointer(buf, totalBufOffset + blockOffsets[curType]), 0, blockLengths[curType]);
390 }
391 }
392 }
393
394 void createIndices(void* buf, size_t bufOffset, void* indices, size_t bufModOffset, int elements) const {
395 int totalIndexOffset = computeActiveElements(bufModOffset); // indices are lineralized and counted up in the loop
396
397 for(int i = 0; i < elements; ++i) {
398 int totalBufOffset = computeBufOffset(i + bufOffset);
399
400 for(int curType = 0; curType < nTypes; ++curType) {
401 types[curType]->createIndices(computeBufferPointer(buf, totalBufOffset + blockOffsets[curType]), 0, indices, totalIndexOffset, blockLengths[curType]);
402 totalIndexOffset += types[curType]->computeActiveElements(blockLengths[curType]);
403 }
404 }
405 }
406
407 void getValues(const void* buf, size_t bufOffset, void* primals, size_t bufModOffset, int elements) const {
408 int totalPrimalsOffset = computeActiveElements(bufModOffset); // indices are lineralized and counted up in the loop
409 for(int i = 0; i < elements; ++i) {
410 int totalBufOffset = computeBufOffset(i + bufOffset);
411
412 for(int curType = 0; curType < nTypes; ++curType) {
413 types[curType]->getValues(computeBufferPointer(buf, totalBufOffset + blockOffsets[curType]), 0, primals, totalPrimalsOffset, blockLengths[curType]);
414 totalPrimalsOffset += types[curType]->computeActiveElements(blockLengths[curType]);
415 }
416 }
417 }
418
419 void performReduce(void* buf, void* target, int count, AMPI_Op op, int ranks) const {
420 for(int j = 1; j < ranks; ++j) {
421 int totalBufOffset = computeBufOffset(count * j);
422
423 MPI_Reduce_local(computeBufferPointer(buf, totalBufOffset), buf, count, this->getMpiType(), op.primalFunction);
424 }
425
426 if(0 != ranks) {
427 copy(buf, 0, target, 0, count);
428 }
429 }
430
431 void copy(void* from, size_t fromOffset, void* to, size_t toOffset, int count) const {
432 for(int i = 0; i < count; ++i) {
433 int totalFromOffset = computeBufOffset(i + fromOffset);
434 int totalToOffset = computeBufOffset(i + toOffset);
435
436 for(int curType = 0; curType < nTypes; ++curType) {
437 types[curType]->copy(computeBufferPointer(from, totalFromOffset + blockOffsets[curType]), 0, computeBufferPointer(to, totalToOffset + modifiedBlockOffsets[curType]), 0, blockLengths[curType]);
438 }
439 }
440 }
441
442 void initializeType(void* buf, size_t bufOffset, int elements) const {
443 for(int i = 0; i < elements; ++i) {
444 int totalBufOffset = computeBufOffset(i + bufOffset);
445
446 for(int curType = 0; curType < nTypes; ++curType) {
447 types[curType]->initializeType(computeBufferPointer(buf, totalBufOffset + blockOffsets[curType]), 0, blockLengths[curType]);
448 }
449 }
450 }
451
452 void freeType(void* buf, size_t bufOffset, int elements) const {
453 for(int i = 0; i < elements; ++i) {
454 int totalBufOffset = computeBufOffset(i + bufOffset);
455
456 for(int curType = 0; curType < nTypes; ++curType) {
457 types[curType]->freeType(computeBufferPointer(buf, totalBufOffset + blockOffsets[curType]), 0, blockLengths[curType]);
458 }
459 }
460 }
461
462 void createTypeBuffer(void* &buf, size_t size) const {
463 char* b = (char*)calloc(size, typeExtent);
464 b -= typeOffset;
465 buf = (void*)b;
466
467 initializeType(buf, 0, size);
468 }
469
470 void createModifiedTypeBuffer(void* &buf, size_t size) const {
471 buf = calloc(size, modifiedExtent);
472 }
473
474
475 void deleteTypeBuffer(void* &buf, size_t size) const {
476 freeType(buf, 0, size);
477
478 char* b = (char*)buf;
479 b += typeOffset;
480 free(b);
481 buf = nullptr;
482 }
483
484 void deleteModifiedTypeBuffer(void* &buf) const {
485 free(buf);
486 buf = nullptr;
487 }
488
490 return new MpiStructType(this);
491 }
492 };
493
494 inline int AMPI_Type_create_contiguous(int count, MpiTypeInterface* oldtype, MpiTypeInterface** newtype) {
495 int typeCount = 1;
496 int* array_of_blocklengths = new int [typeCount];
497 MPI_Aint* array_of_displacements = new MPI_Aint [typeCount];
498 MpiTypeInterface** array_of_types = new MpiTypeInterface*[typeCount];
499
500 array_of_blocklengths[0] = count;
501 array_of_displacements[0] = 0;
502 array_of_types[0] = oldtype;
503
504 *newtype = new MpiStructType(typeCount, array_of_blocklengths, array_of_displacements, array_of_types);
505
506 delete [] array_of_blocklengths;
507 delete [] array_of_displacements;
508 delete [] array_of_types;
509
510 return 0;
511 }
512
513 inline int AMPI_Type_vector(int count, int blocklength, int stride, MpiTypeInterface* oldtype, MpiTypeInterface** newtype) {
514 int typeCount = count;
515 int* array_of_blocklengths = new int [typeCount];
516 MPI_Aint* array_of_displacements = new MPI_Aint [typeCount];
517 MpiTypeInterface** array_of_types = new MpiTypeInterface*[typeCount];
518
519#if MEDI_MPI_TARGET < MEDI_MPI_VERSION_2_0
520 MPI_Aint extent;
521 MPI_Type_extent(oldtype->getMpiType(), &extent);
522#else
523 MPI_Aint extent, lb;
524 MPI_Type_get_extent(oldtype->getMpiType(), &lb, &extent);
525#endif
526
527 for(int i = 0; i < count; ++i) {
528 array_of_blocklengths[i] = blocklength;
529 array_of_displacements[i] = stride * extent * i;
530 array_of_types[i] = oldtype;
531 }
532
533 *newtype = new MpiStructType(typeCount, array_of_blocklengths, array_of_displacements, array_of_types);
534
535 delete [] array_of_blocklengths;
536 delete [] array_of_displacements;
537 delete [] array_of_types;
538
539 return 0;
540 }
541
542 inline int AMPI_Type_create_hvector(int count, int blocklength, MPI_Aint stride, MpiTypeInterface* oldtype, MpiTypeInterface** newtype) {
543 int typeCount = count;
544 int* array_of_blocklengths = new int [typeCount];
545 MPI_Aint* array_of_displacements = new MPI_Aint [typeCount];
546 MpiTypeInterface** array_of_types = new MpiTypeInterface*[typeCount];
547
548 for(int i = 0; i < count; ++i) {
549 array_of_blocklengths[i] = blocklength;
550 array_of_displacements[i] = stride * i;
551 array_of_types[i] = oldtype;
552 }
553
554 *newtype = new MpiStructType(typeCount, array_of_blocklengths, array_of_displacements, array_of_types);
555
556 delete [] array_of_blocklengths;
557 delete [] array_of_displacements;
558 delete [] array_of_types;
559
560 return 0;
561 }
562
563#if MEDI_MPI_TARGET < MEDI_MPI_VERSION_4_0
564 inline int AMPI_Type_hvector(int count, int blocklength, MPI_Aint stride, MpiTypeInterface* oldtype, MpiTypeInterface** newtype) {
565 return AMPI_Type_create_hvector(count, blocklength, stride, oldtype, newtype);
566 }
567#endif
568
569 inline int AMPI_Type_indexed(int count, int* array_of_blocklengths, int* array_of_displacements, MpiTypeInterface* oldtype, MpiTypeInterface** newtype) {
570 int typeCount = count;
571 MPI_Aint* array_of_displacements_byte = new MPI_Aint [typeCount];
572 MpiTypeInterface** array_of_types = new MpiTypeInterface*[typeCount];
573
574#if MEDI_MPI_TARGET < MEDI_MPI_VERSION_2_0
575 MPI_Aint extent;
576 MPI_Type_extent(oldtype->getMpiType(), &extent);
577#else
578 MPI_Aint extent, lb;
579 MPI_Type_get_extent(oldtype->getMpiType(), &lb, &extent);
580#endif
581
582 for(int i = 0; i < count; ++i) {
583 array_of_displacements_byte[i] = array_of_displacements[i] * extent * i;
584 array_of_types[i] = oldtype;
585 }
586
587 *newtype = new MpiStructType(typeCount, array_of_blocklengths, array_of_displacements_byte, array_of_types);
588
589 delete [] array_of_displacements;
590 delete [] array_of_types;
591
592 return 0;
593 }
594
595 inline int AMPI_Type_create_hindexed(int count, int* array_of_blocklengths, MPI_Aint* array_of_displacements, MpiTypeInterface* oldtype, MpiTypeInterface** newtype) {
596 int typeCount = count;
597 MpiTypeInterface** array_of_types = new MpiTypeInterface*[typeCount];
598
599 for(int i = 0; i < count; ++i) {
600 array_of_types[i] = oldtype;
601 }
602
603 *newtype = new MpiStructType(typeCount, array_of_blocklengths, array_of_displacements, array_of_types);
604
605 delete [] array_of_types;
606
607 return 0;
608 }
609
610#if MEDI_MPI_TARGET < MEDI_MPI_VERSION_4_0
611 inline int AMPI_Type_hindexed(int count, int* array_of_blocklength, MPI_Aint* array_of_displacements, MpiTypeInterface* oldtype, MpiTypeInterface** newtype) {
612 return AMPI_Type_create_hindexed(count, array_of_blocklength, array_of_displacements, oldtype, newtype);
613 }
614#endif
615
616 inline int AMPI_Type_indexed_block(int count, int blocklength, int* array_of_displacements, MpiTypeInterface* oldtype, MpiTypeInterface** newtype) {
617 int typeCount = count;
618 int* array_of_blocklengths = new int [typeCount];
619 MPI_Aint* array_of_displacements_byte = new MPI_Aint [typeCount];
620 MpiTypeInterface** array_of_types = new MpiTypeInterface*[typeCount];
621
622#if MEDI_MPI_TARGET < MEDI_MPI_VERSION_2_0
623 MPI_Aint extent;
624 MPI_Type_extent(oldtype->getMpiType(), &extent);
625#else
626 MPI_Aint extent, lb;
627 MPI_Type_get_extent(oldtype->getMpiType(), &lb, &extent);
628#endif
629
630 for(int i = 0; i < count; ++i) {
631 array_of_blocklengths[i] = blocklength;
632 array_of_displacements_byte[i] = array_of_displacements[i] * extent * i;
633 array_of_types[i] = oldtype;
634 }
635
636 *newtype = new MpiStructType(typeCount, array_of_blocklengths, array_of_displacements_byte, array_of_types);
637
638 delete [] array_of_blocklengths;
639 delete [] array_of_displacements;
640 delete [] array_of_types;
641
642 return 0;
643 }
644
645 inline int AMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint* array_of_displacements, MpiTypeInterface* oldtype, MpiTypeInterface** newtype) {
646 int typeCount = count;
647 int* array_of_blocklengths = new int [typeCount];
648 MpiTypeInterface** array_of_types = new MpiTypeInterface*[typeCount];
649
650 for(int i = 0; i < count; ++i) {
651 array_of_blocklengths[i] = blocklength;
652 array_of_types[i] = oldtype;
653 }
654
655 *newtype = new MpiStructType(typeCount, array_of_blocklengths, array_of_displacements, array_of_types);
656
657 delete [] array_of_blocklengths;
658 delete [] array_of_types;
659
660 return 0;
661 }
662
663 inline int dimToOrderDim(int dim, const int dimBase, const int dimStep) {
664 return dim * dimStep + dimBase;
665 }
666
667 inline void add_subarray(const int curDim,
668 const int dimBase,
669 const int dimStep,
670 const int ndims,
671 const int* array_of_subsizes,
672 const int* array_of_starts,
673 int& arrayPos,
674 const MPI_Aint* extents,
675 MPI_Aint dimDisplacement,
676 MPI_Aint* array_of_displacements) {
677 int orderDim = dimToOrderDim(curDim, dimBase, dimStep);
678 if(ndims == curDim + 1) {
679 // I am the last so add the location
680 array_of_displacements[arrayPos] = dimDisplacement + extents[orderDim] * array_of_starts[orderDim];
681 arrayPos += 1;
682 } else {
683 // I am not the last so compute the sub array offset
684 for(int pos = 0; pos < array_of_subsizes[orderDim]; ++pos) {
685 MPI_Aint curDimDisplacement = dimDisplacement + (array_of_starts[orderDim] + pos) * extents[orderDim];
686 add_subarray(curDim + 1, dimBase, dimStep, ndims, array_of_subsizes, array_of_starts, arrayPos, extents,
687 curDimDisplacement, array_of_displacements);
688 }
689 }
690 }
691
692 inline int AMPI_Type_create_subarray(int ndims,
693 const int* array_of_sizes,
694 const int* array_of_subsizes,
695 const int* array_of_starts,
696 int order,
697 MpiTypeInterface* oldtype,
698 MpiTypeInterface** newtype) {
699
700 // decide if to loop from 0 to ndim or ndim to zero
701 int dimBase = 0;
702 int dimStep = 0;
703 if(order == MPI_ORDER_FORTRAN) {
704 dimStep = -1;
705 dimBase = ndims - 1;
706
707 } else if(order == MPI_ORDER_C) {
708 dimStep = 1;
709 dimBase = 0;
710
711 } else {
712 MEDI_EXCEPTION("Unknown order enumerator %d.", order);
713 }
714
715 // compute the total extend of all the blocks
716#if MEDI_MPI_TARGET < MEDI_MPI_VERSION_2_0
717 MPI_Aint extent;
718 MPI_Type_extent(oldtype->getMpiType(), &extent);
719#else
720 MPI_Aint extent, lb;
721 MPI_Type_get_extent(oldtype->getMpiType(), &lb, &extent);
722#endif
723 MPI_Aint* extents = new MPI_Aint [ndims];
724 MPI_Aint curExtent = extent;
725 for(int i = ndims - 1; i >= 0; --i) {
726 int orderDim = dimToOrderDim(i, dimBase, dimStep);
727 extents[orderDim] = curExtent;
728 curExtent *= array_of_sizes[orderDim];
729 }
730
731 // compute the total number of types
732 int typeCount = 1;
733 for(int i = 0; i < ndims - 1; ++i) { // last dimension is used as blocklength
734 typeCount *= array_of_subsizes[dimToOrderDim(i, dimBase, dimStep)];
735 }
736
737 int* array_of_blocklengths = new int [typeCount];
738 MPI_Aint* array_of_displacements = new MPI_Aint [typeCount];
739 MpiTypeInterface** array_of_types = new MpiTypeInterface*[typeCount];
740
741 int arrayPos = 0;
742 add_subarray(0, dimBase, dimStep, ndims, array_of_subsizes, array_of_starts, arrayPos, extents, 0,
743 array_of_displacements);
744
745
746 for(int i = 0; i < typeCount; ++i) {
747 array_of_blocklengths[i] = array_of_subsizes[dimToOrderDim(ndims - 1, dimBase, dimStep)];
748 array_of_types[i] = oldtype;
749 }
750
751 *newtype = new MpiStructType(typeCount, array_of_blocklengths, array_of_displacements, array_of_types);
752
753 delete [] array_of_blocklengths;
754 delete [] array_of_displacements;
755 delete [] array_of_types;
756
757 delete [] extents;
758
759 return 0;
760 }
761
762 inline int AMPI_Type_create_resized(MpiTypeInterface* oldtype, MPI_Aint lb, MPI_Aint extent, MpiTypeInterface** newtype) {
763
764 int typeCount = 1;
765 int* array_of_blocklengths = new int [typeCount];
766 MPI_Aint* array_of_displacements = new MPI_Aint [typeCount];
767 MpiTypeInterface** array_of_types = new MpiTypeInterface*[typeCount];
768
769 array_of_blocklengths[0] = 1;
770 array_of_displacements[0] = 0;
771 array_of_types[0] = oldtype;
772
773 MpiStructType* tempType = new MpiStructType(typeCount, array_of_blocklengths, array_of_displacements, array_of_types);
774
775 *newtype = new MpiStructType(tempType, lb, extent);
776
777 delete tempType;
778
779 delete [] array_of_blocklengths;
780 delete [] array_of_displacements;
781 delete [] array_of_types;
782
783 return 0;
784 }
785
786
787 inline int AMPI_Type_create_struct(int count, MEDI_OPTIONAL_CONST int* array_of_blocklengths, MEDI_OPTIONAL_CONST MPI_Aint* array_of_displacements, MpiTypeInterface* const* array_of_types, MpiTypeInterface** newtype) {
788
789 *newtype = new MpiStructType(count, array_of_blocklengths, array_of_displacements, array_of_types);
790
791 return 0;
792 }
793
794#if MEDI_MPI_TARGET < MEDI_MPI_VERSION_4_0
795 inline int AMPI_Type_struct(int count, MEDI_OPTIONAL_CONST int* array_of_blocklengths, MEDI_OPTIONAL_CONST MPI_Aint* array_of_displacements, MpiTypeInterface* const* array_of_types, MpiTypeInterface** newtype) {
796 return AMPI_Type_create_struct(count, array_of_blocklengths, array_of_displacements, array_of_types, newtype);
797 }
798#endif
799
801 MpiTypeInterface* datatype = *d;
802
803 if(datatype->isModifiedBufferRequired()) {
804 MPI_Datatype modType = datatype->getModifiedMpiType();
805 MPI_Type_commit(&modType);
806 }
807
808 MPI_Datatype type = datatype->getMpiType();
809 return MPI_Type_commit(&type);
810 }
811
813 MpiTypeInterface* datatype = *d;
814
815 delete datatype;
816
817 *d = nullptr;
818 return 0;
819 }
820
821 inline int AMPI_Type_dup(MpiTypeInterface* oldtype, MpiTypeInterface** newtype) {
822
823 *newtype = oldtype->clone();
824
825 return 0;
826 }
827}
The interface for the AD tool that is accessed by MeDiPack.
Definition adToolInterface.h:42
Handling for costum MPI_Datatypes crated by the user.
Definition constructedDatatypes.hpp:48
void PrimalType
Definition constructedDatatypes.hpp:71
void initializeType(void *buf, size_t bufOffset, int elements) const
Initialize the types in the buffer.
Definition constructedDatatypes.hpp:442
void registerValue(void *buf, size_t bufOffset, void *indices, void *oldPrimals, size_t bufModOffset, int elements) const
Register all the AD values on the new machine.
Definition constructedDatatypes.hpp:371
MpiStructType(const MpiStructType *other, size_t offset, size_t extent)
Definition constructedDatatypes.hpp:123
MpiStructType(int count, const int *array_of_blocklengths, const MPI_Aint *array_of_displacements, MpiTypeInterface *const *array_of_types)
Definition constructedDatatypes.hpp:151
void freeType(void *buf, size_t bufOffset, int elements) const
Destroy the types in the buffer.
Definition constructedDatatypes.hpp:452
void copyFromModifiedBuffer(void *buf, size_t bufOffset, const void *bufMod, size_t bufModOffset, int elements) const
Copy all data from the modified buffer and perform the special handling for the AD type.
Definition constructedDatatypes.hpp:346
void AdjointType
Definition constructedDatatypes.hpp:70
void * computeBufferPointer(void *buf, size_t offset) const
Definition constructedDatatypes.hpp:318
MpiStructType(const MpiStructType *other)
Definition constructedDatatypes.hpp:99
void copyIntoModifiedBuffer(const void *buf, size_t bufOffset, void *bufMod, size_t bufModOffset, int elements) const
Copy all data into the modified buffer and perform the special handling for the AD type.
Definition constructedDatatypes.hpp:335
const void * computeBufferPointer(const void *buf, size_t offset) const
Definition constructedDatatypes.hpp:314
void deleteTypeBuffer(void *&buf, size_t size) const
Delete the temporary buffer.
Definition constructedDatatypes.hpp:475
void IndexType
Definition constructedDatatypes.hpp:72
MpiStructType * clone() const
Creates a clone of the mpi type also calling MPI_Type_dub.
Definition constructedDatatypes.hpp:489
bool isModifiedBufferRequired() const
Tell the functions if the underlying AD tool requires new send/recv buffers or if the original buffer...
Definition constructedDatatypes.hpp:322
int computeBufOffset(size_t element) const
Definition constructedDatatypes.hpp:306
void copy(void *from, size_t fromOffset, void *to, size_t toOffset, int count) const
Copy the elements of one buffer into the other.
Definition constructedDatatypes.hpp:431
void createTypeBuffer(void *&buf, size_t size) const
Create a temporary buffer of the type that this interface represents.
Definition constructedDatatypes.hpp:462
int computeActiveElements(const int count) const
Get the number of active elements that are contained in count versions of the type.
Definition constructedDatatypes.hpp:326
void ModifiedType
Definition constructedDatatypes.hpp:69
void performReduce(void *buf, void *target, int count, AMPI_Op op, int ranks) const
Perform a local reduce operation.
Definition constructedDatatypes.hpp:419
int computeModOffset(size_t element) const
Definition constructedDatatypes.hpp:310
void Type
Definition constructedDatatypes.hpp:68
void getValues(const void *buf, size_t bufOffset, void *primals, size_t bufModOffset, int elements) const
Get the primal values from the AD types.
Definition constructedDatatypes.hpp:407
void createIndices(void *buf, size_t bufOffset, void *indices, size_t bufModOffset, int elements) const
Create indices for a reciving buffer if necessary.
Definition constructedDatatypes.hpp:394
void createModifiedTypeBuffer(void *&buf, size_t size) const
Create a temporary buffer of the modified type that this interface represents.
Definition constructedDatatypes.hpp:470
void getIndices(const void *buf, size_t bufOffset, void *indices, size_t bufModOffset, int elements) const
Get all the AD identifiers from all AD types in the buffer.
Definition constructedDatatypes.hpp:358
const ADToolInterface & getADTool() const
Get the AD tool that handled the AD specifics.
Definition constructedDatatypes.hpp:331
void deleteModifiedTypeBuffer(void *&buf) const
Delete the temporary buffer for the modified types.
Definition constructedDatatypes.hpp:484
~MpiStructType()
Definition constructedDatatypes.hpp:284
void clearIndices(void *buf, size_t bufOffset, int elements) const
Clear the AD types in the buffer such that they can be overwritten.
Definition constructedDatatypes.hpp:384
Wrapper interface for MPI types in communications.
Definition typeInterface.hpp:63
virtual void getIndices(const void *buf, size_t bufOffset, void *indices, size_t bufModOffset, int elements) const =0
Get all the AD identifiers from all AD types in the buffer.
virtual void copyFromModifiedBuffer(void *buf, size_t bufOffset, const void *bufMod, size_t bufModOffset, int elements) const =0
Copy all data from the modified buffer and perform the special handling for the AD type.
virtual void freeType(Type *buf, size_t bufOffset, int elements) const =0
Destroy the types in the buffer.
virtual void copy(void *from, size_t fromOffset, void *to, size_t toOffset, int count) const =0
Copy the elements of one buffer into the other.
virtual void getValues(const void *buf, size_t bufOffset, void *primals, size_t bufModOffset, int elements) const =0
Get the primal values from the AD types.
virtual bool isModifiedBufferRequired() const =0
Tell the functions if the underlying AD tool requires new send/recv buffers or if the original buffer...
MPI_Datatype getMpiType() const
Return the MPI type for the data that this interface represents.
Definition typeInterface.hpp:127
MPI_Datatype getModifiedMpiType() const
Return the MPI type for the modified data.
Definition typeInterface.hpp:135
virtual const ADToolInterface & getADTool() const =0
Get the AD tool that handled the AD specifics.
virtual MpiTypeInterface * clone() const =0
Creates a clone of the mpi type also calling MPI_Type_dub.
virtual void registerValue(void *buf, size_t bufOffset, void *indices, void *oldPrimals, size_t bufModOffset, int elements) const =0
Register all the AD values on the new machine.
void setMpiTypes(MPI_Datatype mpiType, MPI_Datatype modifiedMpiType)
Helper method for extending classes to set the types after the initial construction.
Definition typeInterface.hpp:116
virtual void initializeType(Type *buf, size_t bufOffset, int elements) const =0
Initialize the types in the buffer.
virtual void createIndices(void *buf, size_t bufOffset, void *indices, size_t bufModOffset, int elements) const =0
Create indices for a reciving buffer if necessary.
virtual int computeActiveElements(const int count) const =0
Get the number of active elements that are contained in count versions of the type.
virtual void copyIntoModifiedBuffer(const void *buf, size_t bufOffset, void *bufMod, size_t bufModOffset, int elements) const =0
Copy all data into the modified buffer and perform the special handling for the AD type.
virtual void clearIndices(void *buf, size_t bufOffset, int elements) const =0
Clear the AD types in the buffer such that they can be overwritten.
#define MEDI_EXCEPTION(...)
Generates an exception.
Definition exceptions.hpp:44
MpiTypes * mpiTypes
Definition helloWorld.cpp:42
#define MEDI_OPTIONAL_CONST
Definition macros.h:70
#define mediAssert(x)
The assert function for MeDiPack it can be enabled with the preprocessor macro MEDI_EnableAssert=true...
Definition macros.h:91
Global namespace for MeDiPack - Message Differentiation Package.
Definition adjointInterface.hpp:37
int AMPI_Type_create_hindexed(int count, int *array_of_blocklengths, MPI_Aint *array_of_displacements, MpiTypeInterface *oldtype, MpiTypeInterface **newtype)
Definition constructedDatatypes.hpp:595
int AMPI_Type_free(MpiTypeInterface **d)
Definition constructedDatatypes.hpp:812
int AMPI_Type_create_hvector(int count, int blocklength, MPI_Aint stride, MpiTypeInterface *oldtype, MpiTypeInterface **newtype)
Definition constructedDatatypes.hpp:542
int AMPI_Type_dup(MpiTypeInterface *oldtype, MpiTypeInterface **newtype)
Definition constructedDatatypes.hpp:821
int AMPI_Type_create_subarray(int ndims, const int *array_of_sizes, const int *array_of_subsizes, const int *array_of_starts, int order, MpiTypeInterface *oldtype, MpiTypeInterface **newtype)
Definition constructedDatatypes.hpp:692
int AMPI_Type_create_resized(MpiTypeInterface *oldtype, MPI_Aint lb, MPI_Aint extent, MpiTypeInterface **newtype)
Definition constructedDatatypes.hpp:762
int AMPI_Type_commit(MpiTypeInterface **d)
Definition constructedDatatypes.hpp:800
int dimToOrderDim(int dim, const int dimBase, const int dimStep)
Definition constructedDatatypes.hpp:663
void add_subarray(const int curDim, const int dimBase, const int dimStep, const int ndims, const int *array_of_subsizes, const int *array_of_starts, int &arrayPos, const MPI_Aint *extents, MPI_Aint dimDisplacement, MPI_Aint *array_of_displacements)
Definition constructedDatatypes.hpp:667
int AMPI_Type_create_contiguous(int count, MpiTypeInterface *oldtype, MpiTypeInterface **newtype)
Definition constructedDatatypes.hpp:494
int AMPI_Type_create_struct(int count, const int *array_of_blocklengths, const MPI_Aint *array_of_displacements, MpiTypeInterface *const *array_of_types, MpiTypeInterface **newtype)
Definition constructedDatatypes.hpp:787
int AMPI_Type_create_hindexed_block(int count, int blocklength, MPI_Aint *array_of_displacements, MpiTypeInterface *oldtype, MpiTypeInterface **newtype)
Definition constructedDatatypes.hpp:645
int AMPI_Type_indexed_block(int count, int blocklength, int *array_of_displacements, MpiTypeInterface *oldtype, MpiTypeInterface **newtype)
Definition constructedDatatypes.hpp:616
int AMPI_Type_indexed(int count, int *array_of_blocklengths, int *array_of_displacements, MpiTypeInterface *oldtype, MpiTypeInterface **newtype)
Definition constructedDatatypes.hpp:569
int AMPI_Type_vector(int count, int blocklength, int stride, MpiTypeInterface *oldtype, MpiTypeInterface **newtype)
Definition constructedDatatypes.hpp:513
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