File: variablesShapes_hl.cpp

package info (click to toggle)
adios2 2.10.2%2Bdfsg1-3
  • links: PTS, VCS
  • area: main
  • in suites: sid
  • 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 (194 lines) | stat: -rw-r--r-- 6,237 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
187
188
189
190
191
192
193
194
/*
 * Distributed under the OSI-approved Apache License, Version 2.0.  See
 * accompanying file Copyright.txt for details.
 *
 * variablesShapes_hl.cpp : adios2 high-level API example to write and read
 *                   supported Variables shapes using stepping (streaming) mode
 *
 *  Created on: Nov 14, 2019
 *      Author: William F Godoy godoywf@ornl.gov
 */

#include <cstddef>   //std::size_t
#include <iostream>  // std::cout
#include <limits>    // std::numeric_limits
#include <numeric>   //std::iota
#include <stdexcept> //std::exception

#include <adios2.h>
#if ADIOS2_USE_MPI
#include <mpi.h>
#endif

void writer(const std::size_t nx, const std::size_t nsteps, const int rank, const int size)
{
    auto lf_compute = [](const std::size_t step, const std::size_t nx,
                         const int rank) -> std::vector<float> {
        const float value = static_cast<float>(step + rank * nx);
        std::vector<float> array(nx);
        std::iota(array.begin(), array.end(), value);
        return array;
    };

    // You can define variables according to:
    // type <T>: string, uint8_t, int8_t8, ... , float, double
    // shape: global : value or array
    //        local  : value (return a global array), array

    // global shape -> this is the physical dimension across MPI processes
    const adios2::Dims shape = {static_cast<std::size_t>(size * nx)};

    // local start for rank offset -> this is the local origin for the rank
    // domain
    const adios2::Dims start = {static_cast<std::size_t>(rank * nx)};

    // local count -> this is the local size from the local start for the rank
    // domain
    const adios2::Dims count = {nx};

    // adios2::Dims is an alias to std::vector<std::size_t>
    // helps remember the inputs to adios2 functions DefineVariable (write) and
    // SetSelection (read) make sure you always pass std::size_t types

#if ADIOS2_USE_MPI
    adios2::fstream out("variables-shapes_hl.bp", adios2::fstream::out, MPI_COMM_WORLD);
#else
    adios2::fstream out("variables-shapes_hl.bp", adios2::fstream::out);
#endif

    for (size_t step = 0; step < nsteps; ++step)
    {
        // this part mimics the compute portion in an application
        const std::vector<float> array = lf_compute(step, nx, rank);

        // ADIOS2 I/O portion

        // minimize global and local values footprint, by only one rank writing
        // the variables
        if (rank == 0)
        {
            // Global value changing over steps
            out.write("Step", static_cast<uint64_t>(step));

            if (step == 0)
            {
                // Constant Global value
                out.write("GlobalValueString", std::string("ADIOS2 Basics Variable Example"));

                // Constant Local value
                out.write("LocalValueInt32", static_cast<int32_t>(rank), adios2::LocalValue);
            }
        }

        // for this example all ranks write a global and a local array
        out.write("GlobalArray", array.data(), shape, start, count);
        out.write("LocalArray", array.data(), {}, {}, count);

        out.end_step();
    }
    out.close();
}

void reader(const int rank, const int size)
{
    auto lf_ArrayToString = [](const std::vector<float> &array) -> std::string {
        std::string contents = "{ ";
        for (const float value : array)
        {
            contents += std::to_string(static_cast<int>(value)) + " ";
        }
        contents += "}";
        return contents;
    };

// all ranks opening the bp file have access to the entire metadata
#if ADIOS2_USE_MPI
    adios2::fstream in("variables-shapes_hl.bp", adios2::fstream::in, MPI_COMM_WORLD);
#else
    adios2::fstream in("variables-shapes_hl.bp", adios2::fstream::in);
#endif

    // reading in streaming mode, supported by all engines
    // similar to std::getline in std::fstream
    adios2::fstep inStep;
    while (adios2::getstep(in, inStep))
    {
        const std::size_t currentStep = inStep.current_step();

        const std::vector<uint64_t> steps = inStep.read<uint64_t>("Step");
        if (!steps.empty() && rank == 0)
        {
            std::cout << "Found Step " << steps.front() << " in currentStep " << currentStep
                      << "\n";
        }

        const std::vector<std::string> globalValueString =
            inStep.read<std::string>("GlobalValueString");
        if (!globalValueString.empty() && rank == 0)
        {
            std::cout << "Found GlobalValueString " << globalValueString.front()
                      << " in currentStep " << currentStep << "\n";
        }

        const std::vector<int32_t> ranks = inStep.read<int32_t>("Ranks");
        if (!ranks.empty() && rank == 0)
        {
            std::cout << "Found rank " << ranks.front() << " in currentStep " << currentStep
                      << "\n";
        }

        const std::vector<float> globalArray = inStep.read<float>("GlobalArray");
        if (!globalArray.empty() && rank == 0)
        {
            std::cout << "Found globalArray " << lf_ArrayToString(globalArray) + " in currentStep "
                      << currentStep << "\n";
        }

        // default reads block 0
        const std::vector<float> localArray = inStep.read<float>("LocalArray");
        if (!localArray.empty() && rank == 0)
        {
            std::cout << "Found localArray " << lf_ArrayToString(localArray) + " in currentStep "
                      << currentStep << "\n";
        }
        // indicate end of adios2 operations for this step
        in.end_step();
    }
    in.close();
}

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

#if ADIOS2_USE_MPI
    MPI_Comm_rank(MPI_COMM_WORLD, &rank);
    MPI_Comm_size(MPI_COMM_WORLD, &size);
#endif

    try
    {
        constexpr std::size_t nx = 10;
        constexpr std::size_t nsteps = 3;

        writer(nx, nsteps, rank, size);
        reader(rank, size);
    }
    catch (std::exception &e)
    {
        std::cout << "ERROR: ADIOS2 exception: " << e.what() << "\n";
#if ADIOS2_USE_MPI
        MPI_Abort(MPI_COMM_WORLD, -1);
#endif
    }

#if ADIOS2_USE_MPI
    MPI_Finalize();
#endif

    return 0;
}