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