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 GenData(myChars, i, start, count, shape);
73 GenData(myUChars, i, start, count, shape);
74 GenData(myShorts, i, start, count, shape);
75 GenData(myUShorts, i, start, count, shape);
76 GenData(myInts, i, start, count, shape);
77 GenData(myUInts, i, start, count, shape);
78 GenData(myFloats, i, start, count, shape);
79 GenData(myDoubles, i, start, count, shape);
80 GenData(myComplexes, i, start, count, shape);
81 GenData(myDComplexes, i, start, count, shape);
82 engine.Put(bpChars, myChars.data(), adios2::Mode::Sync);
83 engine.Put(bpUChars, myUChars.data(), adios2::Mode::Sync);
84 engine.Put(bpShorts, myShorts.data(), adios2::Mode::Sync);
85 engine.Put(bpUShorts, myUShorts.data(), adios2::Mode::Sync);
86 engine.Put(bpInts, myInts.data(), adios2::Mode::Sync);
87 engine.Put(bpUInts, myUInts.data(), adios2::Mode::Sync);
88 engine.Put(bpFloats, myFloats.data(), adios2::Mode::Sync);
89 engine.Put(bpDoubles, myDoubles.data(), adios2::Mode::Sync);
90 engine.Put(bpComplexes, myComplexes.data(), adios2::Mode::Sync);
91 engine.Put(bpDComplexes, myDComplexes.data(), adios2::Mode::Sync);
92 engine.Put(scalarInt, i);
93 std::string s = "sample string sample string sample string";
94 engine.Put(stringVar, s);
95 engine.EndStep();
96 }
97 engine.Close();
98 }
99
Reader(const Dims & shape,const Dims & start,const Dims & count,const size_t steps,const adios2::Params & engineParams,const std::string & name)100 void Reader(const Dims &shape, const Dims &start, const Dims &count,
101 const size_t steps, const adios2::Params &engineParams,
102 const std::string &name)
103 {
104 adios2::ADIOS adios(mpiComm);
105 adios2::IO dataManIO = adios.DeclareIO("Test");
106 dataManIO.SetEngine("ssc");
107 dataManIO.SetParameters(engineParams);
108 adios2::Engine engine = dataManIO.Open(name, adios2::Mode::Read);
109
110 size_t datasize =
111 std::accumulate(count.begin(), count.end(), static_cast<size_t>(1),
112 std::multiplies<size_t>());
113 std::vector<char> myChars(datasize);
114 std::vector<unsigned char> myUChars(datasize);
115 std::vector<short> myShorts(datasize);
116 std::vector<unsigned short> myUShorts(datasize);
117 std::vector<int> myInts(datasize);
118 std::vector<unsigned int> myUInts(datasize);
119 std::vector<float> myFloats(datasize);
120 std::vector<double> myDoubles(datasize);
121 std::vector<std::complex<float>> myComplexes(datasize);
122 std::vector<std::complex<double>> myDComplexes(datasize);
123
124 engine.LockReaderSelections();
125
126 while (true)
127 {
128 adios2::StepStatus status = engine.BeginStep(StepMode::Read, 5);
129 if (status == adios2::StepStatus::OK)
130 {
131 auto scalarInt = dataManIO.InquireVariable<int>("scalarInt");
132 auto blocksInfo =
133 engine.BlocksInfo(scalarInt, engine.CurrentStep());
134
135 for (const auto &bi : blocksInfo)
136 {
137 ASSERT_EQ(bi.IsValue, true);
138 ASSERT_EQ(bi.Value, engine.CurrentStep());
139 ASSERT_EQ(scalarInt.Min(), engine.CurrentStep());
140 ASSERT_EQ(scalarInt.Max(), engine.CurrentStep());
141 }
142
143 const auto &vars = dataManIO.AvailableVariables();
144 ASSERT_EQ(vars.size(), 12);
145 size_t currentStep = engine.CurrentStep();
146 adios2::Variable<char> bpChars =
147 dataManIO.InquireVariable<char>("bpChars");
148 adios2::Variable<unsigned char> bpUChars =
149 dataManIO.InquireVariable<unsigned char>("bpUChars");
150 adios2::Variable<short> bpShorts =
151 dataManIO.InquireVariable<short>("bpShorts");
152 adios2::Variable<unsigned short> bpUShorts =
153 dataManIO.InquireVariable<unsigned short>("bpUShorts");
154 adios2::Variable<int> bpInts =
155 dataManIO.InquireVariable<int>("bpInts");
156 adios2::Variable<unsigned int> bpUInts =
157 dataManIO.InquireVariable<unsigned int>("bpUInts");
158 adios2::Variable<float> bpFloats =
159 dataManIO.InquireVariable<float>("bpFloats");
160 adios2::Variable<double> bpDoubles =
161 dataManIO.InquireVariable<double>("bpDoubles");
162 adios2::Variable<std::complex<float>> bpComplexes =
163 dataManIO.InquireVariable<std::complex<float>>("bpComplexes");
164 adios2::Variable<std::complex<double>> bpDComplexes =
165 dataManIO.InquireVariable<std::complex<double>>("bpDComplexes");
166 adios2::Variable<std::string> stringVar =
167 dataManIO.InquireVariable<std::string>("stringVar");
168
169 bpChars.SetSelection({start, count});
170 bpUChars.SetSelection({start, count});
171 bpShorts.SetSelection({start, count});
172 bpUShorts.SetSelection({start, count});
173 bpInts.SetSelection({start, count});
174 bpUInts.SetSelection({start, count});
175 bpFloats.SetSelection({start, count});
176 bpDoubles.SetSelection({start, count});
177 bpComplexes.SetSelection({start, count});
178 bpDComplexes.SetSelection({start, count});
179
180 engine.Get(bpChars, myChars.data(), adios2::Mode::Sync);
181 engine.Get(bpUChars, myUChars.data(), adios2::Mode::Sync);
182 engine.Get(bpShorts, myShorts.data(), adios2::Mode::Sync);
183 engine.Get(bpUShorts, myUShorts.data(), adios2::Mode::Sync);
184 engine.Get(bpInts, myInts.data(), adios2::Mode::Sync);
185 engine.Get(bpUInts, myUInts.data(), adios2::Mode::Sync);
186
187 VerifyData(myChars.data(), currentStep, start, count, shape,
188 mpiRank);
189 VerifyData(myUChars.data(), currentStep, start, count, shape,
190 mpiRank);
191 VerifyData(myShorts.data(), currentStep, start, count, shape,
192 mpiRank);
193 VerifyData(myUShorts.data(), currentStep, start, count, shape,
194 mpiRank);
195 VerifyData(myInts.data(), currentStep, start, count, shape,
196 mpiRank);
197 VerifyData(myUInts.data(), currentStep, start, count, shape,
198 mpiRank);
199
200 engine.Get(bpFloats, myFloats.data(), adios2::Mode::Deferred);
201 engine.Get(bpDoubles, myDoubles.data(), adios2::Mode::Deferred);
202 engine.Get(bpComplexes, myComplexes.data(), adios2::Mode::Deferred);
203 engine.Get(bpDComplexes, myDComplexes.data(),
204 adios2::Mode::Deferred);
205 engine.PerformGets();
206
207 VerifyData(myFloats.data(), currentStep, start, count, shape,
208 mpiRank);
209 VerifyData(myDoubles.data(), currentStep, start, count, shape,
210 mpiRank);
211 VerifyData(myComplexes.data(), currentStep, start, count, shape,
212 mpiRank);
213 VerifyData(myDComplexes.data(), currentStep, start, count, shape,
214 mpiRank);
215
216 std::string s;
217 engine.Get(stringVar, s);
218 engine.PerformGets();
219 ASSERT_EQ(s, "sample string sample string sample string");
220 ASSERT_EQ(stringVar.Min(),
221 "sample string sample string sample string");
222 ASSERT_EQ(stringVar.Max(),
223 "sample string sample string sample string");
224
225 int i;
226 engine.Get(scalarInt, &i);
227 engine.PerformGets();
228 ASSERT_EQ(i, currentStep);
229 engine.EndStep();
230 }
231 else if (status == adios2::StepStatus::EndOfStream)
232 {
233 std::cout << "[Rank " + std::to_string(mpiRank) +
234 "] SscTest reader end of stream!"
235 << std::endl;
236 break;
237 }
238 }
239 auto attInt = dataManIO.InquireAttribute<int>("AttInt");
240 std::cout << "[Rank " + std::to_string(mpiRank) + "] Attribute received "
241 << attInt.Data()[0] << ", expected 110" << std::endl;
242 ASSERT_EQ(110, attInt.Data()[0]);
243 ASSERT_NE(111, attInt.Data()[0]);
244 engine.Close();
245 }
246
TEST_F(SscEngineTest,TestSscBase)247 TEST_F(SscEngineTest, TestSscBase)
248 {
249 std::string filename = "TestSscBase";
250 adios2::Params engineParams = {};
251
252 int worldRank, worldSize;
253 MPI_Comm_rank(MPI_COMM_WORLD, &worldRank);
254 MPI_Comm_size(MPI_COMM_WORLD, &worldSize);
255 int mpiGroup = worldRank / (worldSize / 2);
256 MPI_Comm_split(MPI_COMM_WORLD, mpiGroup, worldRank, &mpiComm);
257
258 MPI_Comm_rank(mpiComm, &mpiRank);
259 MPI_Comm_size(mpiComm, &mpiSize);
260
261 Dims shape = {10, (size_t)mpiSize * 2};
262 Dims start = {2, (size_t)mpiRank * 2};
263 Dims count = {5, 2};
264 size_t steps = 100;
265
266 if (mpiGroup == 0)
267 {
268 Writer(shape, start, count, steps, engineParams, filename);
269 }
270
271 std::this_thread::sleep_for(std::chrono::milliseconds(1));
272
273 if (mpiGroup == 1)
274 {
275 Reader(shape, start, count, steps, engineParams, filename);
276 }
277
278 MPI_Barrier(MPI_COMM_WORLD);
279 }
280
main(int argc,char ** argv)281 int main(int argc, char **argv)
282 {
283 MPI_Init(&argc, &argv);
284 int worldRank, worldSize;
285 MPI_Comm_rank(MPI_COMM_WORLD, &worldRank);
286 MPI_Comm_size(MPI_COMM_WORLD, &worldSize);
287 ::testing::InitGoogleTest(&argc, argv);
288 int result = RUN_ALL_TESTS();
289
290 MPI_Finalize();
291 return result;
292 }
293