1 /*
2 * Distributed under the OSI-approved Apache License, Version 2.0. See
3 * accompanying file Copyright.txt for details.
4 */
5
6 #include "TestSscCommon.h"
7 #include <adios2.h>
8 #include <gtest/gtest.h>
9 #include <mpi.h>
10 #include <numeric>
11 #include <thread>
12
13 using namespace adios2;
14 int mpiRank = 0;
15 int mpiSize = 1;
16 MPI_Comm mpiComm;
17
18 class SscEngineTest : public ::testing::Test
19 {
20 public:
21 SscEngineTest() = default;
22 };
23
Writer(const Dims & shape,const Dims & start,const Dims & count,const size_t steps,const adios2::Params & engineParams,const std::string & name)24 void Writer(const Dims &shape, const Dims &start, const Dims &count,
25 const size_t steps, const adios2::Params &engineParams,
26 const std::string &name)
27 {
28 size_t datasize =
29 std::accumulate(count.begin(), count.end(), static_cast<size_t>(1),
30 std::multiplies<size_t>());
31 adios2::ADIOS adios(mpiComm);
32 adios2::IO dataManIO = adios.DeclareIO("WAN");
33 dataManIO.SetEngine("ssc");
34 dataManIO.SetParameters(engineParams);
35 std::vector<char> myChars(datasize);
36 std::vector<unsigned char> myUChars(datasize);
37 std::vector<short> myShorts(datasize);
38 std::vector<unsigned short> myUShorts(datasize);
39 std::vector<int> myInts(datasize);
40 std::vector<unsigned int> myUInts(datasize);
41 std::vector<float> myFloats(datasize);
42 std::vector<double> myDoubles(datasize);
43 std::vector<std::complex<float>> myComplexes(datasize);
44 std::vector<std::complex<double>> myDComplexes(datasize);
45 auto bpChars =
46 dataManIO.DefineVariable<char>("bpChars", shape, start, count);
47 auto bpUChars = dataManIO.DefineVariable<unsigned char>("bpUChars", shape,
48 start, count);
49 auto bpShorts =
50 dataManIO.DefineVariable<short>("bpShorts", shape, start, count);
51 auto bpUShorts = dataManIO.DefineVariable<unsigned short>(
52 "bpUShorts", shape, start, count);
53 auto bpInts = dataManIO.DefineVariable<int>("bpInts", shape, start, count);
54 auto bpUInts =
55 dataManIO.DefineVariable<unsigned int>("bpUInts", shape, start, count);
56 auto bpFloats =
57 dataManIO.DefineVariable<float>("bpFloats", shape, start, count);
58 auto bpDoubles =
59 dataManIO.DefineVariable<double>("bpDoubles", shape, start, count);
60 auto bpComplexes = dataManIO.DefineVariable<std::complex<float>>(
61 "bpComplexes", shape, start, count);
62 auto bpDComplexes = dataManIO.DefineVariable<std::complex<double>>(
63 "bpDComplexes", shape, start, count);
64 auto scalarInt = dataManIO.DefineVariable<int>("scalarInt");
65 auto stringVar = dataManIO.DefineVariable<std::string>("stringVar");
66 dataManIO.DefineAttribute<int>("AttInt", 110);
67 adios2::Engine engine = dataManIO.Open(name, adios2::Mode::Write);
68 // engine.LockWriterDefinitions();
69 for (int i = 0; i < steps; ++i)
70 {
71 engine.BeginStep();
72
73 Dims vstart = {0, (size_t)mpiRank * 2};
74 // Dims vcount = {(size_t)i+1, 2};
75 // Dims vshape = {(size_t)i+1, (size_t)mpiSize * 2};
76 Dims vcount = {10, 2};
77 Dims vshape = {10, (size_t)mpiSize * 2};
78
79 bpChars.SetShape(vshape);
80 bpChars.SetSelection({vstart, vcount});
81 bpUChars.SetShape(vshape);
82 bpUChars.SetSelection({vstart, vcount});
83 bpShorts.SetShape(vshape);
84 bpShorts.SetSelection({vstart, vcount});
85 bpUShorts.SetShape(vshape);
86 bpUShorts.SetSelection({vstart, vcount});
87 bpInts.SetShape(vshape);
88 bpInts.SetSelection({vstart, vcount});
89 bpUInts.SetShape(vshape);
90 bpUInts.SetSelection({vstart, vcount});
91 bpFloats.SetShape(vshape);
92 bpFloats.SetSelection({vstart, vcount});
93 bpDoubles.SetShape(vshape);
94 bpDoubles.SetSelection({vstart, vcount});
95 bpComplexes.SetShape(vshape);
96 bpComplexes.SetSelection({vstart, vcount});
97 bpDComplexes.SetShape(vshape);
98 bpDComplexes.SetSelection({vstart, vcount});
99
100 GenData(myChars, i, vstart, vcount, vshape);
101 GenData(myUChars, i, vstart, vcount, vshape);
102 GenData(myShorts, i, vstart, vcount, vshape);
103 GenData(myUShorts, i, vstart, vcount, vshape);
104 GenData(myInts, i, vstart, vcount, vshape);
105 GenData(myUInts, i, vstart, vcount, vshape);
106 GenData(myFloats, i, vstart, vcount, vshape);
107 GenData(myDoubles, i, vstart, vcount, vshape);
108 GenData(myComplexes, i, vstart, vcount, vshape);
109 GenData(myDComplexes, i, vstart, vcount, vshape);
110
111 engine.Put(bpChars, myChars.data(), adios2::Mode::Sync);
112 engine.Put(bpUChars, myUChars.data(), adios2::Mode::Sync);
113 engine.Put(bpShorts, myShorts.data(), adios2::Mode::Sync);
114 engine.Put(bpUShorts, myUShorts.data(), adios2::Mode::Sync);
115 engine.Put(bpInts, myInts.data(), adios2::Mode::Sync);
116 engine.Put(bpUInts, myUInts.data(), adios2::Mode::Sync);
117 engine.Put(bpFloats, myFloats.data(), adios2::Mode::Sync);
118 engine.Put(bpDoubles, myDoubles.data(), adios2::Mode::Sync);
119 engine.Put(bpComplexes, myComplexes.data(), adios2::Mode::Sync);
120 engine.Put(bpDComplexes, myDComplexes.data(), adios2::Mode::Sync);
121 engine.Put(scalarInt, i);
122 std::string s = "sample string sample string sample string";
123 engine.Put(stringVar, s);
124 engine.EndStep();
125 }
126 engine.Close();
127 }
128
Reader(const Dims & shape,const Dims & start,const Dims & count,const size_t steps,const adios2::Params & engineParams,const std::string & name)129 void Reader(const Dims &shape, const Dims &start, const Dims &count,
130 const size_t steps, const adios2::Params &engineParams,
131 const std::string &name)
132 {
133 adios2::ADIOS adios(mpiComm);
134 adios2::IO dataManIO = adios.DeclareIO("Test");
135 dataManIO.SetEngine("ssc");
136 dataManIO.SetParameters(engineParams);
137 adios2::Engine engine = dataManIO.Open(name, adios2::Mode::Read);
138 // engine.LockReaderSelections();
139
140 size_t datasize =
141 std::accumulate(count.begin(), count.end(), static_cast<size_t>(1),
142 std::multiplies<size_t>());
143 std::vector<char> myChars(datasize);
144 std::vector<unsigned char> myUChars(datasize);
145 std::vector<short> myShorts(datasize);
146 std::vector<unsigned short> myUShorts(datasize);
147 std::vector<int> myInts(datasize);
148 std::vector<unsigned int> myUInts(datasize);
149 std::vector<float> myFloats(datasize);
150 std::vector<double> myDoubles(datasize);
151 std::vector<std::complex<float>> myComplexes(datasize);
152 std::vector<std::complex<double>> myDComplexes(datasize);
153
154 while (true)
155 {
156 adios2::StepStatus status = engine.BeginStep(StepMode::Read, 5);
157 if (status == adios2::StepStatus::OK)
158 {
159 auto scalarInt = dataManIO.InquireVariable<int>("scalarInt");
160 auto blocksInfo =
161 engine.BlocksInfo(scalarInt, engine.CurrentStep());
162
163 for (const auto &bi : blocksInfo)
164 {
165 ASSERT_EQ(bi.IsValue, true);
166 ASSERT_EQ(bi.Value, engine.CurrentStep());
167 ASSERT_EQ(scalarInt.Min(), engine.CurrentStep());
168 ASSERT_EQ(scalarInt.Max(), engine.CurrentStep());
169 }
170
171 const auto &vars = dataManIO.AvailableVariables();
172 ASSERT_EQ(vars.size(), 12);
173 size_t currentStep = engine.CurrentStep();
174 adios2::Variable<char> bpChars =
175 dataManIO.InquireVariable<char>("bpChars");
176 adios2::Variable<unsigned char> bpUChars =
177 dataManIO.InquireVariable<unsigned char>("bpUChars");
178 adios2::Variable<short> bpShorts =
179 dataManIO.InquireVariable<short>("bpShorts");
180 adios2::Variable<unsigned short> bpUShorts =
181 dataManIO.InquireVariable<unsigned short>("bpUShorts");
182 adios2::Variable<int> bpInts =
183 dataManIO.InquireVariable<int>("bpInts");
184 adios2::Variable<unsigned int> bpUInts =
185 dataManIO.InquireVariable<unsigned int>("bpUInts");
186 adios2::Variable<float> bpFloats =
187 dataManIO.InquireVariable<float>("bpFloats");
188 adios2::Variable<double> bpDoubles =
189 dataManIO.InquireVariable<double>("bpDoubles");
190 adios2::Variable<std::complex<float>> bpComplexes =
191 dataManIO.InquireVariable<std::complex<float>>("bpComplexes");
192 adios2::Variable<std::complex<double>> bpDComplexes =
193 dataManIO.InquireVariable<std::complex<double>>("bpDComplexes");
194 adios2::Variable<std::string> stringVar =
195 dataManIO.InquireVariable<std::string>("stringVar");
196
197 auto vshape = bpChars.Shape();
198 myChars.resize(std::accumulate(vshape.begin(), vshape.end(),
199 static_cast<size_t>(1),
200 std::multiplies<size_t>()));
201 myUChars.resize(std::accumulate(vshape.begin(), vshape.end(),
202 static_cast<size_t>(1),
203 std::multiplies<size_t>()));
204 myShorts.resize(std::accumulate(vshape.begin(), vshape.end(),
205 static_cast<size_t>(1),
206 std::multiplies<size_t>()));
207 myUShorts.resize(std::accumulate(vshape.begin(), vshape.end(),
208 static_cast<size_t>(1),
209 std::multiplies<size_t>()));
210 myInts.resize(std::accumulate(vshape.begin(), vshape.end(),
211 static_cast<size_t>(1),
212 std::multiplies<size_t>()));
213 myUInts.resize(std::accumulate(vshape.begin(), vshape.end(),
214 static_cast<size_t>(1),
215 std::multiplies<size_t>()));
216 myFloats.resize(std::accumulate(vshape.begin(), vshape.end(),
217 static_cast<size_t>(1),
218 std::multiplies<size_t>()));
219 myDoubles.resize(std::accumulate(vshape.begin(), vshape.end(),
220 static_cast<size_t>(1),
221 std::multiplies<size_t>()));
222 myComplexes.resize(std::accumulate(vshape.begin(), vshape.end(),
223 static_cast<size_t>(1),
224 std::multiplies<size_t>()));
225 myDComplexes.resize(std::accumulate(vshape.begin(), vshape.end(),
226 static_cast<size_t>(1),
227 std::multiplies<size_t>()));
228
229 engine.Get(bpChars, myChars.data(), adios2::Mode::Sync);
230 engine.Get(bpUChars, myUChars.data(), adios2::Mode::Sync);
231 engine.Get(bpShorts, myShorts.data(), adios2::Mode::Sync);
232 engine.Get(bpUShorts, myUShorts.data(), adios2::Mode::Sync);
233 engine.Get(bpInts, myInts.data(), adios2::Mode::Sync);
234 engine.Get(bpUInts, myUInts.data(), adios2::Mode::Sync);
235 engine.Get(bpFloats, myFloats.data(), adios2::Mode::Sync);
236 engine.Get(bpDoubles, myDoubles.data(), adios2::Mode::Sync);
237 engine.Get(bpComplexes, myComplexes.data(), adios2::Mode::Sync);
238 engine.Get(bpDComplexes, myDComplexes.data(), adios2::Mode::Sync);
239 std::string s;
240 engine.Get(stringVar, s, adios2::Mode::Sync);
241 ASSERT_EQ(s, "sample string sample string sample string");
242 ASSERT_EQ(stringVar.Min(),
243 "sample string sample string sample string");
244 ASSERT_EQ(stringVar.Max(),
245 "sample string sample string sample string");
246
247 int i;
248 engine.Get(scalarInt, &i, adios2::Mode::Sync);
249 ASSERT_EQ(i, currentStep);
250
251 VerifyData(myChars.data(), currentStep, {0, 0}, vshape, vshape,
252 mpiRank);
253 VerifyData(myUChars.data(), currentStep, {0, 0}, vshape, vshape,
254 mpiRank);
255 VerifyData(myShorts.data(), currentStep, {0, 0}, vshape, vshape,
256 mpiRank);
257 VerifyData(myUShorts.data(), currentStep, {0, 0}, vshape, vshape,
258 mpiRank);
259 VerifyData(myInts.data(), currentStep, {0, 0}, vshape, vshape,
260 mpiRank);
261 VerifyData(myUInts.data(), currentStep, {0, 0}, vshape, vshape,
262 mpiRank);
263 VerifyData(myFloats.data(), currentStep, {0, 0}, vshape, vshape,
264 mpiRank);
265 VerifyData(myDoubles.data(), currentStep, {0, 0}, vshape, vshape,
266 mpiRank);
267 VerifyData(myComplexes.data(), currentStep, {0, 0}, vshape, vshape,
268 mpiRank);
269 VerifyData(myDComplexes.data(), currentStep, {0, 0}, vshape, vshape,
270 mpiRank);
271 engine.EndStep();
272 }
273 else if (status == adios2::StepStatus::EndOfStream)
274 {
275 std::cout << "[Rank " + std::to_string(mpiRank) +
276 "] SscTest reader end of stream!"
277 << std::endl;
278 break;
279 }
280 }
281 auto attInt = dataManIO.InquireAttribute<int>("AttInt");
282 std::cout << "[Rank " + std::to_string(mpiRank) + "] Attribute received "
283 << attInt.Data()[0] << ", expected 110" << std::endl;
284 ASSERT_EQ(110, attInt.Data()[0]);
285 ASSERT_NE(111, attInt.Data()[0]);
286 engine.Close();
287 }
288
TEST_F(SscEngineTest,TestSscVaryingSteps)289 TEST_F(SscEngineTest, TestSscVaryingSteps)
290 {
291 std::string filename = "TestSscVaryingSteps";
292 adios2::Params engineParams = {};
293
294 int worldRank, worldSize;
295 MPI_Comm_rank(MPI_COMM_WORLD, &worldRank);
296 MPI_Comm_size(MPI_COMM_WORLD, &worldSize);
297 int mpiGroup = worldRank / (worldSize / 2);
298 MPI_Comm_split(MPI_COMM_WORLD, mpiGroup, worldRank, &mpiComm);
299
300 MPI_Comm_rank(mpiComm, &mpiRank);
301 MPI_Comm_size(mpiComm, &mpiSize);
302
303 Dims shape = {1, (size_t)mpiSize * 2};
304 Dims start = {0, (size_t)mpiRank * 2};
305 Dims count = {1, 2};
306 size_t steps = 100;
307
308 if (mpiGroup == 0)
309 {
310 Writer(shape, start, count, steps, engineParams, filename);
311 }
312
313 std::this_thread::sleep_for(std::chrono::milliseconds(1));
314
315 if (mpiGroup == 1)
316 {
317 Reader(shape, start, count, steps, engineParams, filename);
318 }
319
320 MPI_Barrier(MPI_COMM_WORLD);
321 }
322
main(int argc,char ** argv)323 int main(int argc, char **argv)
324 {
325 MPI_Init(&argc, &argv);
326 int worldRank, worldSize;
327 MPI_Comm_rank(MPI_COMM_WORLD, &worldRank);
328 MPI_Comm_size(MPI_COMM_WORLD, &worldSize);
329 ::testing::InitGoogleTest(&argc, argv);
330 int result = RUN_ALL_TESTS();
331
332 MPI_Finalize();
333 return result;
334 }
335