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