mirror of
https://github.com/OPM/opm-simulators.git
synced 2024-12-22 15:33:29 -06:00
added: simple HDF5File class
bare minimum to enable reading and writing byte arrays
This commit is contained in:
parent
b598c27605
commit
d73d52e162
@ -172,6 +172,9 @@ if(MPI_FOUND)
|
||||
opm/simulators/utils/ParallelSerialization.cpp
|
||||
opm/simulators/utils/SetupZoltanParams.cpp)
|
||||
endif()
|
||||
if(HDF5_FOUND)
|
||||
list(APPEND MAIN_SOURCE_FILES opm/simulators/utils/HDF5File.cpp)
|
||||
endif()
|
||||
|
||||
# originally generated with the command:
|
||||
# find tests -name '*.cpp' -a ! -wholename '*/not-unit/*' -printf '\t%p\n' | sort
|
||||
@ -218,6 +221,9 @@ endif()
|
||||
if(ROCALUTION_FOUND)
|
||||
list(APPEND TEST_SOURCE_FILES tests/test_rocalutionSolver.cpp)
|
||||
endif()
|
||||
if(HDF5_FOUND)
|
||||
list(APPEND TEST_SOURCE_FILES tests/test_HDF5File.cpp)
|
||||
endif()
|
||||
|
||||
list (APPEND TEST_DATA_FILES
|
||||
tests/equil_base.DATA
|
||||
@ -430,6 +436,12 @@ list (APPEND PUBLIC_HEADER_FILES
|
||||
opm/simulators/wells/WGState.hpp
|
||||
)
|
||||
|
||||
if(HDF5_FOUND)
|
||||
list(APPEND PUBLIC_HEADER_FILES
|
||||
opm/simulators/utils/HDF5File.hpp
|
||||
)
|
||||
endif()
|
||||
|
||||
list (APPEND EXAMPLE_SOURCE_FILES
|
||||
examples/printvfp.cpp
|
||||
)
|
||||
|
@ -48,6 +48,7 @@ set (opm-simulators_DEPS
|
||||
"opm-grid REQUIRED"
|
||||
"opm-models REQUIRED"
|
||||
"Damaris 1.7"
|
||||
"HDF5"
|
||||
)
|
||||
|
||||
find_package_deps(opm-simulators)
|
||||
|
126
opm/simulators/utils/HDF5File.cpp
Normal file
126
opm/simulators/utils/HDF5File.cpp
Normal file
@ -0,0 +1,126 @@
|
||||
/*
|
||||
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.
|
||||
*/
|
||||
|
||||
#include <config.h>
|
||||
#include <opm/simulators/utils/HDF5File.hpp>
|
||||
|
||||
#include <filesystem>
|
||||
#include <stdexcept>
|
||||
|
||||
namespace {
|
||||
|
||||
bool groupExists(hid_t parent, const std::string& path)
|
||||
{
|
||||
// turn off errors to avoid cout spew
|
||||
H5E_BEGIN_TRY {
|
||||
#if H5_VERS_MINOR > 8
|
||||
return H5Lexists(parent, path.c_str(), H5P_DEFAULT) == 1;
|
||||
#else
|
||||
return H5Gget_objinfo(static_cast<hid_t>(parent), path.c_str(), 0, nullptr) == 0;
|
||||
#endif
|
||||
} H5E_END_TRY;
|
||||
return false;
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
namespace Opm {
|
||||
|
||||
HDF5File::HDF5File(const std::string& fileName, OpenMode mode)
|
||||
{
|
||||
bool exists = std::filesystem::exists(fileName);
|
||||
if (mode == OpenMode::OVERWRITE ||
|
||||
(mode == OpenMode::APPEND && !exists)) {
|
||||
m_file = H5Fcreate(fileName.c_str(),
|
||||
H5F_ACC_TRUNC,
|
||||
H5P_DEFAULT, H5P_DEFAULT);
|
||||
} else {
|
||||
m_file = H5Fopen(fileName.c_str(),
|
||||
mode == OpenMode::READ ? H5F_ACC_RDONLY : H5F_ACC_RDWR,
|
||||
H5P_DEFAULT);
|
||||
}
|
||||
if (m_file == H5I_INVALID_HID) {
|
||||
throw std::runtime_error(std::string("HDF5File: Failed to ") +
|
||||
( mode == OpenMode::OVERWRITE ||
|
||||
(mode == OpenMode::APPEND && !exists) ? "create" : "open") +
|
||||
fileName);
|
||||
}
|
||||
}
|
||||
|
||||
HDF5File::~HDF5File()
|
||||
{
|
||||
if (m_file != H5I_INVALID_HID) {
|
||||
H5Fclose(m_file);
|
||||
}
|
||||
}
|
||||
|
||||
void HDF5File::write(const std::string& group,
|
||||
const std::string& dset,
|
||||
const std::vector<char>& buffer)
|
||||
{
|
||||
hid_t grp;
|
||||
if (groupExists(m_file, group)) {
|
||||
grp = H5Gopen2(m_file, group.c_str(), H5P_DEFAULT);
|
||||
} else {
|
||||
grp = H5Gcreate2(m_file, group.c_str(), 0, H5P_DEFAULT, H5P_DEFAULT);
|
||||
}
|
||||
|
||||
hsize_t size = buffer.size();
|
||||
hsize_t start = 0;
|
||||
|
||||
hid_t space = H5Screate_simple(1, &size, nullptr);
|
||||
hid_t dataset_id = H5Dcreate2(grp, dset.c_str(), H5T_NATIVE_CHAR, space,
|
||||
H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
|
||||
if (dataset_id == H5I_INVALID_HID) {
|
||||
throw std::runtime_error("HDF5File: Trying to write already existing dataset " + group + '/' + dset);
|
||||
}
|
||||
|
||||
if (size > 0) {
|
||||
hid_t filespace = H5Dget_space(dataset_id);
|
||||
hsize_t stride = 1;
|
||||
H5Sselect_hyperslab(filespace, H5S_SELECT_SET, &start, &stride, &size, nullptr);
|
||||
hid_t memspace = H5Screate_simple(1, &size, nullptr);
|
||||
H5Dwrite(dataset_id, H5T_NATIVE_CHAR, memspace, filespace, H5P_DEFAULT, buffer.data());
|
||||
H5Sclose(memspace);
|
||||
H5Sclose(filespace);
|
||||
}
|
||||
H5Dclose(dataset_id);
|
||||
H5Sclose(space);
|
||||
H5Gclose(grp);
|
||||
}
|
||||
|
||||
void HDF5File::read(const std::string& group,
|
||||
const std::string& dset,
|
||||
std::vector<char>& buffer)
|
||||
{
|
||||
hid_t dataset_id = H5Dopen2(m_file, (group + "/"+ dset).c_str(), H5P_DEFAULT);
|
||||
if (dataset_id == H5I_INVALID_HID) {
|
||||
throw std::runtime_error("HDF5File: Trying to read non-existing dataset " + group + '/' + dset);
|
||||
}
|
||||
|
||||
hid_t space = H5Dget_space(dataset_id);
|
||||
hsize_t size = H5Sget_simple_extent_npoints(space);
|
||||
buffer.resize(size);
|
||||
H5Dread(dataset_id, H5T_NATIVE_CHAR, H5S_ALL, H5S_ALL, H5P_DEFAULT, buffer.data());
|
||||
H5Dclose(dataset_id);
|
||||
}
|
||||
|
||||
}
|
73
opm/simulators/utils/HDF5File.hpp
Normal file
73
opm/simulators/utils/HDF5File.hpp
Normal file
@ -0,0 +1,73 @@
|
||||
/*
|
||||
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 HDF5_FILE_HPP
|
||||
#define HDF5_FILE_HPP
|
||||
|
||||
#include <hdf5.h>
|
||||
|
||||
#include <string>
|
||||
#include <vector>
|
||||
|
||||
namespace Opm {
|
||||
|
||||
//! \brief Class handling simple output to HDF5.
|
||||
class HDF5File {
|
||||
public:
|
||||
//! \brief Enumeration of file opening modes.
|
||||
enum class OpenMode {
|
||||
APPEND, //!< Append to an existing file (creates new if not)
|
||||
OVERWRITE, //!< Overwrite and write to an existing file
|
||||
READ //!< Open existing file for reading
|
||||
};
|
||||
|
||||
//! \brief Opens HDF5 file for I/O.
|
||||
//! \param fileName Name of file to open
|
||||
//! \param mode Open mode for file
|
||||
HDF5File(const std::string& fileName, OpenMode mode);
|
||||
|
||||
//! \brief Destructor clears up any opened files.
|
||||
~HDF5File();
|
||||
|
||||
//! \brief Write a char buffer to a specified location in file.
|
||||
//! \param group Group ("directory") to write data to
|
||||
//! \param dset Data set ("file") to write data to
|
||||
//! \param buffer Data to write
|
||||
//! \details Throws exception on failure
|
||||
void write(const std::string& group,
|
||||
const std::string& dset,
|
||||
const std::vector<char>& buffer);
|
||||
|
||||
//! \brief Read a char buffer from a specified location in file.
|
||||
//! \param group Group ("directory") to read data from
|
||||
//! \param dset Data set ("file") to read data from
|
||||
//! \param buffer Vector to store read data in
|
||||
//! \details Throws exception on failure
|
||||
void read(const std::string& group,
|
||||
const std::string& dset,
|
||||
std::vector<char>& buffer);
|
||||
|
||||
private:
|
||||
hid_t m_file = H5I_INVALID_HID; //!< File handle
|
||||
};
|
||||
|
||||
}
|
||||
|
||||
#endif
|
96
tests/test_HDF5File.cpp
Normal file
96
tests/test_HDF5File.cpp
Normal file
@ -0,0 +1,96 @@
|
||||
/*
|
||||
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 <opm/simulators/utils/HDF5File.hpp>
|
||||
|
||||
#define BOOST_TEST_MODULE HDF5FileTest
|
||||
#include <boost/test/unit_test.hpp>
|
||||
|
||||
#include <filesystem>
|
||||
#include <stdexcept>
|
||||
|
||||
using namespace Opm;
|
||||
|
||||
BOOST_AUTO_TEST_CASE(ReadWrite)
|
||||
{
|
||||
auto path = std::filesystem::temp_directory_path() / Opm::unique_path("hdf5test%%%%%");
|
||||
std::filesystem::create_directory(path);
|
||||
auto rwpath = (path / "rw.hdf5").string();
|
||||
const std::vector<char> test_data{1,2,3,4,5,6,8,9};
|
||||
{
|
||||
Opm::HDF5File out_file(rwpath, Opm::HDF5File::OpenMode::OVERWRITE);
|
||||
BOOST_CHECK_NO_THROW(out_file.write("/test_data", "d1", test_data));
|
||||
}
|
||||
{
|
||||
Opm::HDF5File in_file(rwpath, Opm::HDF5File::OpenMode::READ);
|
||||
std::vector<char> data;
|
||||
BOOST_CHECK_NO_THROW(in_file.read("/test_data", "d1", data));
|
||||
BOOST_CHECK_EQUAL_COLLECTIONS(data.begin(), data.end(),
|
||||
test_data.begin(), test_data.end());
|
||||
}
|
||||
std::filesystem::remove(rwpath);
|
||||
std::filesystem::remove(path);
|
||||
}
|
||||
|
||||
BOOST_AUTO_TEST_CASE(ThrowOpenNonexistent)
|
||||
{
|
||||
BOOST_CHECK_THROW(Opm::HDF5File out_file("no_such_file.hdf5", Opm::HDF5File::OpenMode::READ), std::runtime_error);
|
||||
}
|
||||
|
||||
BOOST_AUTO_TEST_CASE(ReadNonExistentDset)
|
||||
{
|
||||
auto path = std::filesystem::temp_directory_path() / Opm::unique_path("hdf5test%%%%%");
|
||||
std::filesystem::create_directory(path);
|
||||
auto rwpath = (path / "existent_dset.hdf5").string();
|
||||
const std::vector<char> test_data{1,2,3,4,5,6,8,9};
|
||||
{
|
||||
Opm::HDF5File out_file(rwpath, Opm::HDF5File::OpenMode::OVERWRITE);
|
||||
BOOST_CHECK_NO_THROW(out_file.write("/test_data", "d1", test_data));
|
||||
}
|
||||
{
|
||||
Opm::HDF5File in_file(rwpath, Opm::HDF5File::OpenMode::READ);
|
||||
std::vector<char> data;
|
||||
BOOST_CHECK_NO_THROW(in_file.read("/test_data", "d1", data));
|
||||
BOOST_CHECK_EQUAL_COLLECTIONS(data.begin(), data.end(),
|
||||
test_data.begin(), test_data.end());
|
||||
BOOST_CHECK_THROW(in_file.read("/test_data", "d2", data), std::runtime_error);
|
||||
}
|
||||
std::filesystem::remove(rwpath);
|
||||
std::filesystem::remove(path);
|
||||
}
|
||||
|
||||
BOOST_AUTO_TEST_CASE(WriteExistentDset)
|
||||
{
|
||||
auto path = std::filesystem::temp_directory_path() / Opm::unique_path("hdf5test%%%%%");
|
||||
std::filesystem::create_directory(path);
|
||||
auto rwpath = (path / "existent_dset.hdf5").string();
|
||||
const std::vector<char> test_data{1,2,3,4,5,6,8,9};
|
||||
{
|
||||
Opm::HDF5File out_file(rwpath, Opm::HDF5File::OpenMode::OVERWRITE);
|
||||
BOOST_CHECK_NO_THROW(out_file.write("/test_data", "d1", test_data));
|
||||
BOOST_CHECK_THROW(out_file.write("/test_data", "d1", test_data), std::runtime_error);
|
||||
}
|
||||
std::filesystem::remove(rwpath);
|
||||
std::filesystem::remove(path);
|
||||
}
|
||||
|
Loading…
Reference in New Issue
Block a user