-
Notifications
You must be signed in to change notification settings - Fork 16
/
variables-shapes.cpp
275 lines (233 loc) · 9 KB
/
variables-shapes.cpp
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
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
/*
* Distributed under the OSI-approved Apache License, Version 2.0. See
* accompanying file Copyright.txt for details.
*
* * variables-shapes.cpp : adios2 low-level API example to write and read
* supported Variables shapes using stepping
* (streaming) mode
*
* Created on: Nov 14, 2019
* Author: William F Godoy [email protected]
*/
#include <cstddef> //std::size_t
#include <iostream> // std::cout
#include <limits> // std::numeric_limits
#include <numeric> //std::iota
#include <stdexcept> //std::exception
#include <adios2.h>
#if ADIOS2_USE_MPI
#include <mpi.h>
#endif
void writer(adios2::ADIOS &adios, const std::size_t nx,
const std::size_t nsteps, const int rank, const int size)
{
auto lf_compute = [](const std::size_t step, const std::size_t nx,
const int rank) -> std::vector<float> {
const float value = static_cast<float>(step + rank * nx);
std::vector<float> array(nx);
std::iota(array.begin(), array.end(), value);
return array;
};
adios2::IO io = adios.DeclareIO("variables-shapes_writer");
// You can define variables according to:
// type <T>: string, uint8_t, int8_t8, ... , float, double
// shape: global : value or array
// local : value (return a global array), array
/********** GLOBAL VALUE *************/
// string variables are always of value type, can't pass dimensions
adios2::Variable<std::string> varGlobalValueString =
io.DefineVariable<std::string>("GlobalValueString");
// global value can change on each step. Example: Step
adios2::Variable<uint64_t> varStep = io.DefineVariable<uint64_t>("Step");
/********** GLOBAL ARRAYS *************/
// For a regular 1D decomposition:
// 0*nx 1*nx 2*nx 3*nx shape (4*nx)
//--------//-------//--------//---------//
// nx nx nx nx
// global shape -> this is the physical dimension across MPI processes
const adios2::Dims shape = {static_cast<std::size_t>(size * nx)};
// local start for rank offset -> this is the local origin for the rank
// domain
const adios2::Dims start = {static_cast<std::size_t>(rank * nx)};
// local count -> this is the local size from the local start for the rank
// domain
const adios2::Dims count = {nx};
// adios2::Dims is an alias to std::vector<std::size_t>
// helps remember the inputs to adios2 functions DefineVariable (write) and
// SetSelection (read) make sure you always pass std::size_t types
// global array with dimensions (shape, start, count)
// last argument indicates to adios2 that dimensions won't change
adios2::Variable<float> varGlobalArray = io.DefineVariable<float>(
"GlobalArray", shape, start, count, adios2::ConstantDims);
/********** LOCAL VALUE **************/
// Independent values per rank at write, but presented as a global array at
// read Example: store current rank, but presented as an array of ranks
adios2::Variable<int32_t> varLocalValueInt32 =
io.DefineVariable<int32_t>("Rank", {adios2::LocalValueDim});
/********** LOCAL ARRAY **************/
// Independent values and dimensions per rank, there is no notion of
// "continuity", each start from 0-origin to count. Example: mesh
// nodes-per-rank
adios2::Variable<float> varLocalArray = io.DefineVariable<float>(
"LocalArray", {}, {}, count, adios2::ConstantDims);
adios2::Engine writer = io.Open("variables-shapes.bp", adios2::Mode::Write);
for (size_t step = 0; step < nsteps; ++step)
{
// this part mimics the compute portion in an application
const std::vector<float> array = lf_compute(step, nx, rank);
// ADIOS2 I/O portion
// BeginStep/EndStep is the streaming API -> supported by all engines
writer.BeginStep();
// minimize global and local values footprint, by only one rank putting
// the variables
if (rank == 0)
{
// Global value changing over steps
writer.Put(varStep, static_cast<uint64_t>(step));
if (step == 0)
{
// Global absolute value
writer.Put(varGlobalValueString,
std::string("ADIOS2 Basics Variable Example"));
// Local absolute value
writer.Put(varLocalValueInt32, static_cast<int32_t>(rank));
}
}
// for this example all ranks put a global and a local array
writer.Put(varGlobalArray, array.data());
writer.Put(varLocalArray, array.data());
writer.EndStep();
}
writer.Close();
}
void reader(adios2::ADIOS &adios, const int rank, const int size)
{
adios2::IO io = adios.DeclareIO("variables-shapes_reader");
// all ranks opening the bp file have access to the entire metadata
adios2::Engine reader = io.Open("variables-shapes.bp", adios2::Mode::Read);
// reading in streaming mode
while (reader.BeginStep() != adios2::StepStatus::EndOfStream)
{
// scope between BeginStep and EndStep is only for the current step
const size_t currentStep = reader.CurrentStep();
// Typical flow: InquireVariable
adios2::Variable<uint64_t> varStep =
io.InquireVariable<uint64_t>("Step");
uint64_t step = std::numeric_limits<uint64_t>::max();
// check Variable existence
if (varStep)
{
if (rank == 0)
{
// variable objects are "printable" reporting Name and Type
std::cout << "Found Global Value " << varStep << " in step "
<< currentStep << "\n";
// output: Found Global Value Variable<uint64_t>(Name: "Step")
// in step 0
}
reader.Get(varStep, step);
}
// GlobalValueString
adios2::Variable<std::string> varGlobalValueString =
io.InquireVariable<std::string>("GlobalValueString");
std::string globalValueString;
// check Variable existence and Get
if (varGlobalValueString)
{
if (rank == 0)
{
std::cout << "Found Global Value " << varGlobalValueString
<< " in step " << currentStep << "\n";
}
reader.Get(varGlobalValueString, globalValueString);
}
// Global Arrays at read from local values at write
adios2::Variable<int32_t> varRanks =
io.InquireVariable<int32_t>("Ranks");
std::vector<int32_t> ranks;
if (varRanks)
{
if (rank == 0)
{
std::cout << "Found Global Array " << varRanks << " in step "
<< currentStep << "\n";
}
// passing a vector convenience: adios2 would resize it
// automatically
reader.Get(varRanks, ranks);
}
// Global Array
adios2::Variable<float> varGlobalArray =
io.InquireVariable<float>("GlobalArray");
std::vector<float> globalArray;
if (varGlobalArray)
{
if (rank == 0)
{
std::cout << "Found GlobalArray " << varGlobalArray
<< " in step " << currentStep << "\n";
}
reader.Get(varGlobalArray, globalArray);
}
// Local Array
adios2::Variable<float> varLocalArray =
io.InquireVariable<float>("LocalArray");
std::vector<float> localArray;
if (varLocalArray)
{
// local arrays require an extra step to select the block of
// interest (0 is default) we only select block 0 in this example
varLocalArray.SetBlockSelection(0);
if (rank == 0)
{
std::cout << "Found LocalArray " << varLocalArray << " in step "
<< currentStep << "\n";
}
reader.Get(varLocalArray, localArray);
}
// since all Get calls are "deferred" all the data would be populated at
// EndStep
reader.EndStep();
// data is available
if (rank == 0)
{
std::cout << "\n";
}
}
reader.Close();
}
int main(int argc, char *argv[])
{
#if ADIOS2_USE_MPI
MPI_Init(&argc, &argv);
#endif
int rank = 0;
int size = 1;
#if ADIOS2_USE_MPI
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
MPI_Comm_size(MPI_COMM_WORLD, &size);
#endif
try
{
#if ADIOS2_USE_MPI
adios2::ADIOS adios(MPI_COMM_WORLD);
#else
adios2::ADIOS adios;
#endif
constexpr std::size_t nx = 10;
constexpr std::size_t nsteps = 3;
writer(adios, nx, nsteps, rank, size);
reader(adios, rank, size);
}
catch (const std::exception &e)
{
std::cout << "ERROR: ADIOS2 exception: " << e.what() << "\n";
#if ADIOS2_USE_MPI
MPI_Abort(MPI_COMM_WORLD, -1);
#endif
}
#if ADIOS2_USE_MPI
MPI_Finalize();
#endif
return 0;
}