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
|
/*
* Distributed under the OSI-approved Apache License, Version 2.0. See
* accompanying file Copyright.txt for details.
*
* hdf5Reader.cpp
*
* Created on: Jan 24, 2018
* Author: Junmin Gu
*/
#include <ios> //std::ios_base::failure
#include <iostream> //std::cout
#include <mpi.h>
#include <stdexcept> //std::invalid_argument std::exception
#include <vector>
#include <adios2.h>
template <class T>
void ReadData(adios2::IO h5IO, adios2::Engine &h5Reader, const std::string &name)
{
adios2::Variable<T> var = h5IO.InquireVariable<T>(name);
if (var)
{
size_t nDims = var.Shape().size();
size_t totalSize = 1;
for (size_t i = 0; i < nDims; i++)
{
totalSize *= var.Shape()[i];
}
std::vector<T> myValues(totalSize);
// myFloats.data is pre-allocated
h5Reader.Get<T>(var, myValues.data(), adios2::Mode::Sync);
// std::cout << "\tValues of "<<name<<": ";
std::cout << "\tPeek Values: ";
if (totalSize < 20)
{ // print all
for (const auto number : myValues)
{
std::cout << number << " ";
}
}
else
{
size_t counter = 0;
for (const auto number : myValues)
{
if ((counter < 5) || (counter > totalSize - 5))
{
std::cout << number << " ";
}
else if (counter == 5)
{
std::cout << " ...... ";
}
counter++;
}
}
std::cout << "\n";
}
}
int main(int argc, char *argv[])
{
int provided;
// MPI_THREAD_MULTIPLE is only required if you enable the SST MPI_DP
MPI_Init_thread(&argc, &argv, MPI_THREAD_MULTIPLE, &provided);
int rank, size;
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
MPI_Comm_size(MPI_COMM_WORLD, &size);
/** Application variable */
const std::size_t Nx = 10;
std::vector<float> myFloats(Nx);
std::vector<int> myInts(Nx);
const char *filename = "myVector.h5";
if (argc > 1)
{
filename = argv[1];
}
std::cout << " Using file: " << filename << std::endl;
try
{
/** ADIOS class factory of IO class objects */
adios2::ADIOS adios(MPI_COMM_WORLD);
/*** IO class object: settings and factory of Settings: Variables,
* Parameters, Transports, and Execution: Engines */
adios2::IO h5IO = adios.DeclareIO("ReadHDF5");
h5IO.SetEngine("HDF5");
/** Engine derived class, spawned to start IO operations */
adios2::Engine h5Reader = h5IO.Open(filename, adios2::Mode::Read);
h5Reader.BeginStep();
const std::map<std::string, adios2::Params> variables = h5IO.AvailableVariables();
for (const auto &variablePair : variables)
{
std::cout << "Name: " << variablePair.first;
std::cout << std::endl;
for (const auto ¶meter : variablePair.second)
{
std::cout << "\t" << parameter.first << ": " << parameter.second << "\n";
if (parameter.second == "double")
{
ReadData<double>(h5IO, h5Reader, variablePair.first);
}
else if (parameter.second == "float")
{
ReadData<float>(h5IO, h5Reader, variablePair.first);
}
else if (parameter.second == "unsigned int")
{
ReadData<unsigned int>(h5IO, h5Reader, variablePair.first);
}
else if (parameter.second == "int")
{
ReadData<int>(h5IO, h5Reader, variablePair.first);
}
//... add more types if needed
}
}
const std::map<std::string, adios2::Params> attributes = h5IO.AvailableAttributes();
for (const auto &attrPair : attributes)
{
std::cout << "AttrName: " << attrPair.first;
std::cout << std::endl;
for (const auto ¶meter : attrPair.second)
{
std::cout << "\t" << parameter.first << ": " << parameter.second << "\n";
if (parameter.second == "double")
{
// ReadData<double>(h5IO, h5Reader, variablePair.first);
}
else if (parameter.second == "float")
{
// ReadData<float>(h5IO, h5Reader, variablePair.first);
}
else if (parameter.second == "unsigned int")
{
// ReadData<unsigned int>(h5IO, h5Reader,
// variablePair.first);
}
else if (parameter.second == "int")
{
// ReadData<int>(h5IO, h5Reader, variablePair.first);
}
//... add more types if needed
}
}
h5Reader.EndStep();
h5Reader.Close();
}
catch (std::invalid_argument &e)
{
std::cout << "Invalid argument exception, STOPPING PROGRAM from rank " << rank << "\n";
std::cout << e.what() << "\n";
}
catch (std::ios_base::failure &e)
{
std::cout << "IO System base failure exception, STOPPING PROGRAM "
"from rank "
<< rank << "\n";
std::cout << e.what() << "\n";
}
catch (std::exception &e)
{
std::cout << "Exception, STOPPING PROGRAM from rank " << rank << "\n";
std::cout << e.what() << "\n";
}
MPI_Finalize();
return 0;
}
|