mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-01-14 04:41:56 -06:00
added: HDF5Serializer
this takes the serialization data and stores it in a hdf5 file, alternatively reads the data and deserializes. will be used for restarting purposes
This commit is contained in:
parent
a90189c78d
commit
4b07f6d010
@ -225,6 +225,7 @@ if(ROCALUTION_FOUND)
|
||||
endif()
|
||||
if(HDF5_FOUND)
|
||||
list(APPEND TEST_SOURCE_FILES tests/test_HDF5File.cpp)
|
||||
list(APPEND TEST_SOURCE_FILES tests/test_HDF5Serializer.cpp)
|
||||
endif()
|
||||
|
||||
list (APPEND TEST_DATA_FILES
|
||||
@ -441,6 +442,7 @@ list (APPEND PUBLIC_HEADER_FILES
|
||||
|
||||
if(HDF5_FOUND)
|
||||
list(APPEND PUBLIC_HEADER_FILES
|
||||
ebos/hdf5serializer.hh
|
||||
opm/simulators/utils/HDF5File.hpp
|
||||
)
|
||||
endif()
|
||||
|
132
ebos/hdf5serializer.hh
Normal file
132
ebos/hdf5serializer.hh
Normal file
@ -0,0 +1,132 @@
|
||||
/*
|
||||
This file is part of the Open Porous Media project (OPM).
|
||||
|
||||
OPM is free software: you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation, either version 2 of the License, or
|
||||
(at your option) any later version.
|
||||
|
||||
OPM is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License
|
||||
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
||||
|
||||
Consult the COPYING file in the top-level source directory of this
|
||||
module for the precise wording of the license and the list of
|
||||
copyright holders.
|
||||
*/
|
||||
#ifndef ECL_HDF5_SERIALIZER_HH
|
||||
#define ECL_HDF5_SERIALIZER_HH
|
||||
|
||||
#include <opm/common/utility/Serializer.hpp>
|
||||
#include <opm/common/utility/MemPacker.hpp>
|
||||
|
||||
#include <opm/simulators/utils/HDF5File.hpp>
|
||||
#include <opm/simulators/utils/moduleVersion.hpp>
|
||||
|
||||
#include <algorithm>
|
||||
#include <cctype>
|
||||
#include <cstdlib>
|
||||
|
||||
namespace Opm {
|
||||
|
||||
//! \brief Class for (de-)serializing using HDF5.
|
||||
class HDF5Serializer : public Serializer<Serialization::MemPacker> {
|
||||
public:
|
||||
HDF5Serializer(const std::string& fileName, HDF5File::OpenMode mode)
|
||||
: Serializer<Serialization::MemPacker>(m_packer_priv)
|
||||
, m_h5file(fileName, mode)
|
||||
{}
|
||||
|
||||
//! \brief Serialize and write data to restart file.
|
||||
//! \tparam T Type of class to write
|
||||
//! \param data Class to write restart data for
|
||||
template<class T>
|
||||
void write(T& data,
|
||||
const std::string& group,
|
||||
const std::string& dset)
|
||||
{
|
||||
try {
|
||||
this->pack(data);
|
||||
} catch (...) {
|
||||
m_packSize = std::numeric_limits<size_t>::max();
|
||||
throw;
|
||||
}
|
||||
|
||||
m_h5file.write(group, dset, m_buffer);
|
||||
}
|
||||
|
||||
//! \brief Writes a header to the file.
|
||||
//! \param simulator_name Name of simulator used
|
||||
//! \param module_version Version of simulator used
|
||||
//! \param time_stamp Build time-stamp for simulator used
|
||||
//! \param case_name Name of case file is associated with
|
||||
//! \param params List of parameter values
|
||||
//! \param num_procs Number of processes used
|
||||
void writeHeader(const std::string& simulator_name,
|
||||
const std::string& module_version,
|
||||
const std::string& time_stamp,
|
||||
const std::string& case_name,
|
||||
const std::string& params,
|
||||
int num_procs)
|
||||
{
|
||||
try {
|
||||
this->pack(simulator_name, module_version, time_stamp,
|
||||
case_name, params, num_procs);
|
||||
} catch (...) {
|
||||
m_packSize = std::numeric_limits<size_t>::max();
|
||||
throw;
|
||||
}
|
||||
m_h5file.write("/", "simulator_info", m_buffer);
|
||||
}
|
||||
|
||||
//! \brief Read data and deserialize from restart file.
|
||||
//! \tparam T Type of class to read
|
||||
//! \param data Class to read restart data for
|
||||
template<class T>
|
||||
void read(T& data,
|
||||
const std::string& group,
|
||||
const std::string& dset)
|
||||
{
|
||||
m_h5file.read(group, dset, m_buffer);
|
||||
this->unpack(data);
|
||||
}
|
||||
|
||||
//! \brief Returns the last report step stored in file.
|
||||
int lastReportStep() const
|
||||
{
|
||||
const auto entries = m_h5file.list("/report_step");
|
||||
int last = -1;
|
||||
for (const auto& entry : entries) {
|
||||
int num = std::atoi(entry.c_str());
|
||||
last = std::max(last, num);
|
||||
}
|
||||
|
||||
return last;
|
||||
}
|
||||
|
||||
//! \brief Returns a list of report steps stored in restart file.
|
||||
std::vector<int> reportSteps() const
|
||||
{
|
||||
const auto entries = m_h5file.list("/report_step");
|
||||
std::vector<int> result(entries.size());
|
||||
std::transform(entries.begin(), entries.end(), result.begin(),
|
||||
[](const std::string& input)
|
||||
{
|
||||
return std::atoi(input.c_str());
|
||||
});
|
||||
std::sort(result.begin(), result.end());
|
||||
return result;
|
||||
}
|
||||
|
||||
private:
|
||||
const Serialization::MemPacker m_packer_priv{}; //!< Packer instance
|
||||
HDF5File m_h5file; //!< HDF5 backend for the serializer
|
||||
};
|
||||
|
||||
}
|
||||
|
||||
#endif
|
83
tests/test_HDF5Serializer.cpp
Normal file
83
tests/test_HDF5Serializer.cpp
Normal file
@ -0,0 +1,83 @@
|
||||
/*
|
||||
Copyright 2021 Equinor.
|
||||
|
||||
This file is part of the Open Porous Media project (OPM).
|
||||
|
||||
OPM is free software: you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation, either version 3 of the License, or
|
||||
(at your option) any later version.
|
||||
|
||||
OPM is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License
|
||||
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
|
||||
#include <config.h>
|
||||
|
||||
#include <opm/common/utility/FileSystem.hpp>
|
||||
|
||||
#include <ebos/hdf5serializer.hh>
|
||||
|
||||
#include <opm/input/eclipse/Schedule/Group/Group.hpp>
|
||||
|
||||
#define BOOST_TEST_MODULE HDF5FileTest
|
||||
#include <boost/test/unit_test.hpp>
|
||||
|
||||
#include <filesystem>
|
||||
#include <stdexcept>
|
||||
|
||||
using namespace Opm;
|
||||
|
||||
BOOST_AUTO_TEST_CASE(Header)
|
||||
{
|
||||
auto path = std::filesystem::temp_directory_path() / Opm::unique_path("hdf5test%%%%%");
|
||||
std::filesystem::create_directory(path);
|
||||
auto rwpath = (path / "rw.hdf5").string();
|
||||
std::array<std::string,5> output{"foo", "bar", "foobar", "bob", "bobbar"};
|
||||
{
|
||||
HDF5Serializer ser(rwpath, HDF5File::OpenMode::OVERWRITE);
|
||||
ser.writeHeader(output[0], output[1], output[2], output[3], output[4], 5);
|
||||
}
|
||||
{
|
||||
HDF5Serializer ser(rwpath, HDF5File::OpenMode::READ);
|
||||
std::tuple<std::array<std::string,5>,int> input;
|
||||
ser.read(input, "/", "simulator_info");
|
||||
const auto& [strings, num_procs] = input;
|
||||
BOOST_CHECK_EQUAL_COLLECTIONS(strings.begin(), strings.end(),
|
||||
output.begin(), output.end());
|
||||
BOOST_CHECK_EQUAL(num_procs, 5);
|
||||
}
|
||||
|
||||
std::filesystem::remove(rwpath);
|
||||
std::filesystem::remove(path);
|
||||
}
|
||||
|
||||
BOOST_AUTO_TEST_CASE(WriteRead)
|
||||
{
|
||||
auto path = std::filesystem::temp_directory_path() / Opm::unique_path("hdf5test%%%%%");
|
||||
std::filesystem::create_directory(path);
|
||||
auto rwpath = (path / "rw.hdf5").string();
|
||||
auto output = Group::serializationTestObject();
|
||||
{
|
||||
HDF5Serializer ser(rwpath, HDF5File::OpenMode::OVERWRITE);
|
||||
ser.write(output, "/report_step/10", "test");
|
||||
}
|
||||
{
|
||||
HDF5Serializer ser(rwpath, HDF5File::OpenMode::READ);
|
||||
Group input;
|
||||
ser.read(input, "/report_step/10", "test");
|
||||
BOOST_CHECK_MESSAGE(input == output, "Deserialized data does not match input");
|
||||
BOOST_CHECK_EQUAL(ser.lastReportStep(), 10);
|
||||
const auto steps = ser.reportSteps();
|
||||
BOOST_CHECK_EQUAL(steps.size(), 1u);
|
||||
BOOST_CHECK_EQUAL(steps[0], 10);
|
||||
}
|
||||
|
||||
std::filesystem::remove(rwpath);
|
||||
std::filesystem::remove(path);
|
||||
}
|
Loading…
Reference in New Issue
Block a user