Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
/*
* Distributed under the OSI-approved Apache License, Version 2.0. See
* accompanying file Copyright.txt for details.
*
* adiosMPIFunctions.tcc : specialization of template functions defined in
* adiosMPIFunctions.h
*
* Created on: Aug 30, 2017
* Author: William F Godoy godoywf@ornl.gov
*/
#ifndef ADIOS2_HELPER_ADIOSMPIFUNCTIONS_TCC_
#define ADIOS2_HELPER_ADIOSMPIFUNCTIONS_TCC_
#include "adiosMPIFunctions.h"
#include <algorithm> //std::foreach
#include <iostream> //TODO: remove
#include <numeric> //std::accumulate
#include "adios2/ADIOSMPI.h"
#include "adios2/ADIOSTypes.h"
#include "adios2/helper/adiosType.h"
namespace adios2
{
// BroadcastValue specializations
template <>
size_t BroadcastValue(const size_t &input, MPI_Comm mpiComm,
const int rankSource)
{
int rank;
MPI_Comm_rank(mpiComm, &rank);
size_t output = 0;
if (rank == rankSource)
{
output = input;
}
MPI_Bcast(&output, 1, ADIOS2_MPI_SIZE_T, rankSource, mpiComm);
return output;
}
template <>
std::string BroadcastValue(const std::string &input, MPI_Comm mpiComm,
const int rankSource)
{
int rank;
MPI_Comm_rank(mpiComm, &rank);
const size_t inputSize = input.size();
size_t length = BroadcastValue(inputSize, mpiComm, rankSource);
std::string output;
if (rank == rankSource)
{
output = input;
}
else
{
output.resize(length);
}
MPI_Bcast(const_cast<char *>(output.data()), length, MPI_CHAR, rankSource,
mpiComm);
return output;
}
// ReduceValue specializations
template <>
unsigned int ReduceValues(const unsigned int source, MPI_Comm mpiComm,
MPI_Op operation, const int rankDestination)
{
const unsigned int sourceLocal = source;
unsigned int reduceValue = 0;
MPI_Reduce(&sourceLocal, &reduceValue, 1, MPI_UNSIGNED, operation,
rankDestination, mpiComm);
return reduceValue;
}
template <>
unsigned long int ReduceValues(const unsigned long int source, MPI_Comm mpiComm,
MPI_Op operation, const int rankDestination)
{
const unsigned long int sourceLocal = source;
unsigned long int reduceValue = 0;
MPI_Reduce(&sourceLocal, &reduceValue, 1, MPI_UNSIGNED_LONG, operation,
rankDestination, mpiComm);
return reduceValue;
}
template <>
unsigned long long int ReduceValues(const unsigned long long int source,
MPI_Comm mpiComm, MPI_Op operation,
const int rankDestination)
{
const unsigned long long int sourceLocal = source;
unsigned long long int reduceValue = 0;
MPI_Reduce(&sourceLocal, &reduceValue, 1, MPI_UNSIGNED_LONG_LONG, operation,
rankDestination, mpiComm);
return reduceValue;
}
// BroadcastVector specializations
template <>
std::vector<char> BroadcastVector(const std::vector<char> &input,
MPI_Comm mpiComm, const int rankSource)
{
// First Broadcast the size, then the contents
size_t inputSize = BroadcastValue(input.size(), mpiComm, rankSource);
int rank;
MPI_Comm_rank(mpiComm, &rank);
std::vector<char> output;
if (rank == rankSource)
{
MPI_Bcast(const_cast<char *>(input.data()), static_cast<int>(inputSize),
MPI_CHAR, rankSource, mpiComm);
return input; // no copy
}
else
{
output.resize(inputSize);
MPI_Bcast(output.data(), static_cast<int>(inputSize), MPI_CHAR,
rankSource, mpiComm);
}
return output;
}
template <>
std::vector<size_t> BroadcastVector(const std::vector<size_t> &input,
MPI_Comm mpiComm, const int rankSource)
{
// First Broadcast the size, then the contents
size_t inputSize = BroadcastValue(input.size(), mpiComm, rankSource);
int rank;
MPI_Comm_rank(mpiComm, &rank);
std::vector<size_t> output;
if (rank == rankSource)
{
MPI_Bcast(const_cast<size_t *>(input.data()),
static_cast<int>(inputSize), ADIOS2_MPI_SIZE_T, rankSource,
mpiComm);
return input; // no copy in rankSource
}
else
{
output.resize(inputSize);
MPI_Bcast(output.data(), static_cast<int>(inputSize), ADIOS2_MPI_SIZE_T,
rankSource, mpiComm);
}
return output;
}
// GatherArrays specializations
template <>
void GatherArrays(const char *source, const size_t sourceCount,
char *destination, MPI_Comm mpiComm,
const int rankDestination)
{
const int countsInt = static_cast<int>(sourceCount);
int result = MPI_Gather(source, countsInt, MPI_CHAR, destination, countsInt,
MPI_CHAR, rankDestination, mpiComm);
if (result != MPI_SUCCESS)
{
throw std::runtime_error("ERROR: in ADIOS2 detected failure in MPI "
"Gather type MPI_CHAR function\n");
}
}
template <>
void GatherArrays(const size_t *source, const size_t sourceCount,
size_t *destination, MPI_Comm mpiComm,
const int rankDestination)
{
const int countsInt = static_cast<int>(sourceCount);
int result =
MPI_Gather(source, countsInt, ADIOS2_MPI_SIZE_T, destination, countsInt,
ADIOS2_MPI_SIZE_T, rankDestination, mpiComm);
if (result != MPI_SUCCESS)
{
throw std::runtime_error("ERROR: in ADIOS2 detected failure in MPI "
"Gather type size_t function\n");
}
}
// GathervArrays specializations
template <>
void GathervArrays(const char *source, const size_t sourceCount,
const size_t *counts, const size_t countsSize,
char *destination, MPI_Comm mpiComm,
const int rankDestination)
{
int result = 0;
int rank;
MPI_Comm_rank(mpiComm, &rank);
std::vector<int> countsInt, displacementsInt;
if (rank == rankDestination)
{
countsInt = NewVectorTypeFromArray<size_t, int>(counts, countsSize);
displacementsInt = GetGathervDisplacements(counts, countsSize);
}
result = MPI_Gatherv(source, static_cast<int>(sourceCount), MPI_CHAR,
destination, countsInt.data(), displacementsInt.data(),
MPI_CHAR, rankDestination, mpiComm);
if (result != MPI_SUCCESS)
{
throw std::runtime_error("ERROR: in ADIOS2 detected failure in MPI "
"Gatherv type MPI_CHAR function\n");
}
}
template <>
void GathervArrays(const size_t *source, const size_t sourceCount,
const size_t *counts, const size_t countsSize,
size_t *destination, MPI_Comm mpiComm,
const int rankDestination)
{
int result = 0;
int rank;
MPI_Comm_rank(mpiComm, &rank);
std::vector<int> countsInt =
NewVectorTypeFromArray<size_t, int>(counts, countsSize);
std::vector<int> displacementsInt =
GetGathervDisplacements(counts, countsSize);
result = MPI_Gatherv(source, sourceCount, ADIOS2_MPI_SIZE_T, destination,
countsInt.data(), displacementsInt.data(),
ADIOS2_MPI_SIZE_T, rankDestination, mpiComm);
if (result != MPI_SUCCESS)
{
throw std::runtime_error("ERROR: in ADIOS2 detected failure in MPI "
"Gather type size_t function\n");
}
}
} // end namespace adios2
#endif /* ADIOS2_HELPER_ADIOSMPIFUNCTIONS_TCC_ */