File: inlineReaderWriter.cpp

package info (click to toggle)
adios2 2.10.2%2Bdfsg1-2
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 33,764 kB
  • sloc: cpp: 175,964; ansic: 160,510; f90: 14,630; yacc: 12,668; python: 7,275; perl: 7,126; sh: 2,825; lisp: 1,106; xml: 1,049; makefile: 579; lex: 557
file content (186 lines) | stat: -rw-r--r-- 6,065 bytes parent folder | download | duplicates (2)
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
/*
 * Distributed under the OSI-approved Apache License, Version 2.0.  See
 * accompanying file Copyright.txt for details.
 *
 * inlineReaderWriter.cpp  example borrowed from bpTimeWriter, using
 * the inline engine. Writes a variable using the Advance function for time
 * aggregation. Time step is saved as an additional (global) single value
 * variable, just for tracking purposes.
 *
 *  Created on: Nov 16, 2018
 *      Author: Aron Helser aron.helser@kitware.com
 */

#include <algorithm> //std::for_each
#include <ios>       //std::ios_base::failure
#include <iostream>  //std::cout
#include <stdexcept> //std::invalid_argument std::exception
#include <vector>

#include <adios2.h>

void DoAnalysis(adios2::IO &inlineIO, adios2::Engine &inlineReader, int rank, unsigned int step)
{
    inlineReader.BeginStep();
    /////////////////////READ
    adios2::Variable<float> inlineFloats000 = inlineIO.InquireVariable<float>("inlineFloats000");

    adios2::Variable<std::string> inlineString =
        inlineIO.InquireVariable<std::string>("inlineString");

    if (inlineFloats000)
    {
        auto blocksInfo = inlineReader.BlocksInfo(inlineFloats000, step);

        std::cout << "Data StepsStart " << inlineFloats000.StepsStart() << " from rank " << rank
                  << ": ";
        for (auto &info : blocksInfo)
        {
            // bp file reader would see all blocks, inline only sees local
            // writer's block(s).
            size_t myBlock = info.BlockID;
            inlineFloats000.SetBlockSelection(myBlock);

            // info passed by reference
            // engine must remember data pointer (or info) to fill it out at
            // PerformGets()
            inlineReader.Get<float>(inlineFloats000, info, adios2::Mode::Deferred);
        }
        inlineReader.PerformGets();

        for (const auto &info : blocksInfo)
        {
            adios2::Dims count = info.Count;
            const float *vectData = info.Data();
            for (size_t i = 0; i < count[0]; ++i)
            {
                float datum = vectData[i];
                std::cout << datum << " ";
            }
            std::cout << "\n";
        }
    }
    else
    {
        std::cout << "Variable inlineFloats000 not found\n";
    }

    if (inlineString && rank == 0)
    {
        std::string myString;
        inlineReader.Get(inlineString, myString, adios2::Mode::Sync);
        std::cout << "inlineString: " << myString << "\n";
    }
    inlineReader.EndStep();
    // all deferred block info are now valid - need data pointers to be
    // valid, filled with data
}

int main(int argc, char *argv[])
{
    int rank = 0, size = 1;
#if ADIOS2_USE_MPI
    int provided;

    // MPI_THREAD_MULTIPLE is only required if you enable the SST MPI_DP
    MPI_Init_thread(&argc, &argv, MPI_THREAD_MULTIPLE, &provided);
    MPI_Comm_rank(MPI_COMM_WORLD, &rank);
    MPI_Comm_size(MPI_COMM_WORLD, &size);
    adios2::ADIOS adios(MPI_COMM_WORLD);
#else
    adios2::ADIOS adios;
#endif

    // Application variable
    std::vector<float> myFloats = {0, 1, 2, 3, 4, 5, 6, 7, 8, 9};
    const std::size_t Nx = myFloats.size();

    try
    {
        // Inline uses single IO for write/read
        adios2::IO inlineIO = adios.DeclareIO("InlineReadWrite");
        /// WRITE
        inlineIO.SetEngine("Inline");
        inlineIO.SetParameter("verbose", "4");

        /** global array: name, { shape (total dimensions) }, { start
         * (local) },
         * { count (local) }, all are constant dimensions */
        const unsigned int variablesSize = 10;
        std::vector<adios2::Variable<float>> inlineFloats(variablesSize);

        adios2::Variable<std::string> inlineString =
            inlineIO.DefineVariable<std::string>("inlineString");

        for (unsigned int v = 0; v < variablesSize; ++v)
        {
            std::string namev("inlineFloats");
            if (v < 10)
            {
                namev += "00";
            }
            else if (v < 100)
            {
                namev += "0";
            }
            namev += std::to_string(v);

            inlineFloats[v] = inlineIO.DefineVariable<float>(namev, {size * Nx}, {rank * Nx}, {Nx},
                                                             adios2::ConstantDims);
        }

        /** global single value variable: name */
        adios2::Variable<unsigned int> inlineTimeStep =
            inlineIO.DefineVariable<unsigned int>("timeStep");

        adios2::Engine inlineWriter = inlineIO.Open("myWriteID", adios2::Mode::Write);

        adios2::Engine inlineReader = inlineIO.Open("myReadID", adios2::Mode::Read);

        for (unsigned int timeStep = 0; timeStep < 3; ++timeStep)
        {
            inlineWriter.BeginStep();
            if (rank == 0) // global single value, only saved by rank 0
            {
                inlineWriter.Put<unsigned int>(inlineTimeStep, timeStep);
            }

            // template type is optional, but recommended
            for (unsigned int v = 0; v < variablesSize; ++v)
            {
                // Note: Put is deferred, so all variables will see v == 9
                // and myFloats[0] == 9, 10, or 11
                myFloats[rank] = static_cast<float>(v + timeStep + rank);
                inlineWriter.Put(inlineFloats[v], myFloats.data());
            }

            const std::string myString("Hello from rank: " + std::to_string(rank) +
                                       " and timestep: " + std::to_string(timeStep));

            if (rank == 0)
            {
                inlineWriter.Put(inlineString, myString);
            }

            inlineWriter.EndStep();

            DoAnalysis(inlineIO, inlineReader, rank, timeStep);
        }
    }
    catch (std::exception const &e)
    {
        std::cout << "Caught exception from rank " << rank << "\n";
        std::cout << e.what() << "\n";
#if ADIOS2_USE_MPI
        return MPI_Abort(MPI_COMM_WORLD, 1);
#else
        return 1;
#endif
    }

#if ADIOS2_USE_MPI
    MPI_Finalize();
#endif

    return 0;
}