diff --git a/src/mesh/ParallelChecker.hpp b/src/mesh/ParallelChecker.hpp index 0c82d8750edfdd9dc13ef85abf1c678293078207..05f90bba9c0b6c5e87418bdf98c2286075af47a7 100644 --- a/src/mesh/ParallelChecker.hpp +++ b/src/mesh/ParallelChecker.hpp @@ -74,7 +74,7 @@ class ParallelChecker }(); auto values_group_id = HDF5::createOrOpenGroup(file_id, "/values"); - auto group_id = HDF5::createOrOpenSubGroup(values_group_id, std::to_string(m_tag)); + auto group_id = HDF5::createOrOpenGroup(values_group_id, std::to_string(m_tag)); HDF5::writeAttribute(group_id, "filename", std::string{source_location.file_name()}); HDF5::writeAttribute(group_id, "function", source_location.function_name()); @@ -86,22 +86,22 @@ class ParallelChecker HDF5::writeAttribute(group_id, "item_type", itemName(item_type)); HDF5::writeAttribute(group_id, "data_type", demangle<DataType>()); - HDF5::write(group_id, name, item_value.arrayView()); + HDF5::writeArray(group_id, name, item_value.arrayView()); switch (i_connectivity->dimension()) { case 1: { const Connectivity<1>& connectivity = dynamic_cast<const Connectivity<1>&>(*i_connectivity); - HDF5::write(group_id, "numbers", connectivity.number<item_type>().arrayView()); + HDF5::writeArray(group_id, "numbers", connectivity.number<item_type>().arrayView()); break; } case 2: { const Connectivity<2>& connectivity = dynamic_cast<const Connectivity<2>&>(*i_connectivity); - HDF5::write(group_id, "numbers", connectivity.number<item_type>().arrayView()); + HDF5::writeArray(group_id, "numbers", connectivity.number<item_type>().arrayView()); break; } case 3: { const Connectivity<3>& connectivity = dynamic_cast<const Connectivity<3>&>(*i_connectivity); - HDF5::write(group_id, "numbers", connectivity.number<item_type>().arrayView()); + HDF5::writeArray(group_id, "numbers", connectivity.number<item_type>().arrayView()); break; } default: { @@ -111,6 +111,8 @@ class ParallelChecker ++m_tag; + HDF5::close(values_group_id); + HDF5::close(group_id); HDF5::close(file_id); std::cout << rang::fg::cyan << " | writing " << rang::fgB::green << "success" << rang::fg::reset << '\n'; @@ -127,7 +129,7 @@ class ParallelChecker auto file_id = HDF5::openFileRO(m_filename); auto values_group_id = HDF5::openGroup(file_id, "/values"); - auto group_id = HDF5::openSubGroup(values_group_id, std::to_string(m_tag)); + auto group_id = HDF5::openGroup(values_group_id, std::to_string(m_tag)); const std::string reference_name = HDF5::readAttribute<std::string>(group_id, "name"); const std::string reference_file_name = HDF5::readAttribute<std::string>(group_id, "filename"); @@ -305,6 +307,8 @@ class ParallelChecker } } + HDF5::close(values_group_id); + HDF5::close(group_id); HDF5::close(file_id); ++m_tag; } diff --git a/src/utils/HDF5.cpp b/src/utils/HDF5.cpp index a0a57d57da2843e414e4c4b5231b5ed4a77fd458..1361ddbe9afd6133fa069b8ce05cfba68d0cc8a8 100644 --- a/src/utils/HDF5.cpp +++ b/src/utils/HDF5.cpp @@ -4,6 +4,7 @@ #include <utils/Exceptions.hpp> +#include <fstream> #include <rang.hpp> #include <sstream> @@ -14,7 +15,7 @@ HDF5::create(const std::string& filename) if (file_id == -1) { std::ostringstream error_msg; - error_msg << "Cannot create file '" << rang::fgB::yellow << filename << rang::fg::reset << "'"; + error_msg << "cannot create file '" << rang::fgB::yellow << filename << rang::fg::reset << "'"; throw NormalError(error_msg.str()); } return file_id; @@ -25,9 +26,9 @@ HDF5::openFileRW(const std::string& filename) { FileId file_id = H5Fopen(filename.c_str(), H5F_ACC_RDWR, H5P_DEFAULT); - if (file_id == -1) { + if (file_id < 0) { std::ostringstream error_msg; - error_msg << "Cannot open file '" << rang::fgB::yellow << filename << rang::fg::reset << "' for writing"; + error_msg << "cannot open file '" << rang::fgB::yellow << filename << rang::fg::reset << "' for writing"; throw NormalError(error_msg.str()); } return file_id; @@ -38,9 +39,9 @@ HDF5::openFileRO(const std::string& filename) { FileId file_id = H5Fopen(filename.c_str(), H5F_ACC_RDONLY, H5P_DEFAULT); - if (file_id == -1) { + if (file_id < 0) { std::ostringstream error_msg; - error_msg << "Cannot read file '" << rang::fgB::yellow << filename << rang::fg::reset << "'"; + error_msg << "cannot open file '" << rang::fgB::yellow << filename << rang::fg::reset << "' for reading"; throw NormalError(error_msg.str()); } return file_id; @@ -61,7 +62,7 @@ HDF5::createOrOpenGroup(const FileId& file_id, const std::string& groupname) } HDF5::GroupId -HDF5::createOrOpenSubGroup(const GroupId& group_id, const std::string& groupname) +HDF5::createOrOpenGroup(const GroupId& group_id, const std::string& groupname) { GroupId sub_group_id = [&] { if (H5Lexists(group_id, groupname.c_str(), H5P_DEFAULT) > 0) { @@ -81,30 +82,30 @@ HDF5::openGroup(const FileId& file_id, const std::string& groupname) return H5Gopen2(file_id, groupname.c_str(), H5P_DEFAULT); } else { std::ostringstream error_msg; - error_msg << "Cannot read group '" << rang::fgB::yellow << groupname << rang::fg::reset << "'"; + error_msg << "cannot open group '" << rang::fgB::yellow << groupname << rang::fg::reset << "'"; throw NormalError(error_msg.str()); } } HDF5::GroupId -HDF5::openSubGroup(const GroupId& group_id, const std::string& groupname) +HDF5::openGroup(const GroupId& group_id, const std::string& groupname) { if (H5Lexists(group_id, groupname.c_str(), H5P_DEFAULT) > 0) { return H5Gopen2(group_id, groupname.c_str(), H5P_DEFAULT); } else { std::ostringstream error_msg; - error_msg << "Cannot read sub-group '" << rang::fgB::yellow << groupname << rang::fg::reset << "'"; + error_msg << "cannot open group '" << rang::fgB::yellow << groupname << rang::fg::reset << "'"; throw NormalError(error_msg.str()); } } template <HDF5::IdType id_type> void -HDF5::_writeDataset(const HDFId<id_type>& id, - const std::string& name, - const DatatypeId& type, - const size_t& array_size, - const void* array) +HDF5::_writeArray(const HDFId<id_type>& id, + const std::string& name, + const DatatypeId& type, + const size_t& array_size, + const void* array) { hsize_t dataspace_size[] = {array_size}; hid_t dataspace = H5Screate_simple(1, dataspace_size, NULL); @@ -112,28 +113,68 @@ HDF5::_writeDataset(const HDFId<id_type>& id, hid_t dataset = H5Dcreate2(id, name.c_str(), type, dataspace, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT); if (dataset == -1) { - throw NormalError("Cannot create dataset"); + throw NormalError("cannot create dataset"); // LCOV_EXCL_LINE } herr_t status = H5Dwrite(dataset, type, H5S_ALL, H5S_ALL, H5P_DEFAULT, array); if (status < 0) { - throw NormalError("Cannot write dataset"); + throw NormalError("cannot write dataset"); // LCOV_EXCL_LINE + } + + H5Dclose(dataset); + H5Sclose(dataspace); +} + +template void HDF5::_writeArray(const FileId& id, + const std::string& name, + const DatatypeId& type, + const size_t& array_size, + const void* array); + +template void HDF5::_writeArray(const GroupId& id, + const std::string& name, + const DatatypeId& type, + const size_t& array_size, + const void* array); + +template <HDF5::IdType id_type> +void +HDF5::_writeTable(const HDFId<id_type>& id, + const std::string& name, + const DatatypeId& type, + const size_t& table_nb_rows, + const size_t& table_nb_columns, + const void* table) +{ + hsize_t dataspace_size[] = {table_nb_rows, table_nb_columns}; + hid_t dataspace = H5Screate_simple(2, dataspace_size, NULL); + + hid_t dataset = H5Dcreate2(id, name.c_str(), type, dataspace, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT); + + if (dataset == -1) { + throw NormalError("cannot create dataset"); // LCOV_EXCL_LINE + } + herr_t status = H5Dwrite(dataset, type, H5S_ALL, H5S_ALL, H5P_DEFAULT, table); + if (status < 0) { + throw NormalError("cannot write dataset"); // LCOV_EXCL_LINE } H5Dclose(dataset); H5Sclose(dataspace); } -template void HDF5::_writeDataset(const FileId& id, - const std::string& name, - const DatatypeId& type, - const size_t& array_size, - const void* array); +template void HDF5::_writeTable(const FileId& id, + const std::string& name, + const DatatypeId& type, + const size_t& table_nb_rows, + const size_t& table_nb_columns, + const void* table); -template void HDF5::_writeDataset(const GroupId& id, - const std::string& name, - const DatatypeId& type, - const size_t& array_size, - const void* array); +template void HDF5::_writeTable(const GroupId& id, + const std::string& name, + const DatatypeId& type, + const size_t& table_nb_rows, + const size_t& table_nb_columns, + const void* table); template <HDF5::IdType id_type> void @@ -141,9 +182,17 @@ HDF5::_writeAttribute(const HDFId<id_type>& id, const std::string& name, const D { hid_t space = H5Screate(H5S_SCALAR); - hid_t attribute = H5Acreate2(id, name.c_str(), type, space, H5P_DEFAULT, H5P_DEFAULT); + if (H5Aexists(id, name.c_str()) > 0) { + std::ostringstream error_msg; + error_msg << "attribute '" << rang::fgB::yellow << name << rang::fg::reset << "' already exists"; + throw NormalError(error_msg.str()); + } + + hid_t attribute_id = H5Acreate2(id, name.c_str(), type, space, H5P_DEFAULT, H5P_DEFAULT); + if (attribute_id < 0) { + } - H5Awrite(attribute, type, value); + H5Awrite(attribute_id, type, value); } template void HDF5::_writeAttribute(const FileId& id, diff --git a/src/utils/HDF5.hpp b/src/utils/HDF5.hpp index c0e03503aa42c743e78c5cf85023c4000df548e8..933d1658d7f6dfe6d48af8610f27049e3110e60e 100644 --- a/src/utils/HDF5.hpp +++ b/src/utils/HDF5.hpp @@ -7,6 +7,7 @@ #include <utils/Array.hpp> #include <utils/Exceptions.hpp> +#include <utils/Table.hpp> #include <hdf5.h> @@ -25,7 +26,7 @@ class HDF5 }; template <IdType type> - class HDFId + class [[nodiscard]] HDFId { private: hid_t m_hdf_id; @@ -70,11 +71,19 @@ class HDF5 private: template <IdType id_type> - static void _writeDataset(const HDFId<id_type>& id, - const std::string& data_name, - const DatatypeId& data_type, - const size_t& array_size, - const void* array); + static void _writeArray(const HDFId<id_type>& id, + const std::string& data_name, + const DatatypeId& data_type, + const size_t& array_size, + const void* array); + + template <IdType id_type> + static void _writeTable(const HDFId<id_type>& id, + const std::string& data_name, + const DatatypeId& data_type, + const size_t& table_nb_rows, + const size_t& table_nb_columns, + const void* table); template <IdType id_type> static void _writeAttribute(const HDFId<id_type>& id, @@ -173,6 +182,16 @@ class HDF5 DatatypeId datatype = H5Tcopy(H5T_NATIVE_LDOUBLE); H5Tset_order(datatype, H5T_ORDER_LE); return datatype; + } else if constexpr (is_tiny_vector_v<std::decay_t<T>>) { + using ContentT = typename T::data_type; + hsize_t dim[] = {DataT::Dimension}; + DatatypeId datatype = H5Tarray_create(_getDataType<ContentT>(), 1, dim); + return datatype; + } else if constexpr (is_tiny_matrix_v<std::decay_t<T>>) { + using ContentT = typename T::data_type; + hsize_t dim[] = {DataT::NumberOfRows, DataT::NumberOfColumns}; + DatatypeId datatype = H5Tarray_create(_getDataType<ContentT>(), 2, dim); + return datatype; } else { // static_assert(std::is_same_v<int, T>, "unexpected data type"); @@ -185,24 +204,30 @@ class HDF5 static FileId openFileRO(const std::string& file_name); static GroupId createOrOpenGroup(const FileId& file_id, const std::string& group_name); - static GroupId createOrOpenSubGroup(const GroupId& group_id, const std::string& group_name); + static GroupId createOrOpenGroup(const GroupId& group_id, const std::string& group_name); static GroupId openGroup(const FileId& file_id, const std::string& group_name); - static GroupId openSubGroup(const GroupId& group_id, const std::string& group_name); + static GroupId openGroup(const GroupId& group_id, const std::string& group_name); template <typename DataT, IdType id_type> static void - write(const HDFId<id_type>& id, const std::string& name, const Array<DataT>& array) + writeArray(const HDFId<id_type>& id, const std::string& name, const Array<DataT>& array) { + if (H5Lexists(id, name.c_str(), H5P_DEFAULT)) { + std::ostringstream error_msg; + error_msg << "dataset '" << name << "' already exists"; + throw NormalError(error_msg.str()); + } + static_assert((id_type == IdType::file) or (id_type == IdType::group)); DatatypeId datatype = _getDataType<DataT>(); H5Tset_order(datatype, H5T_ORDER_LE); const void* array_start = &(array[typename Array<DataT>::index_type(0)]); - _writeDataset(id, name, datatype, array.size(), array_start); + _writeArray(id, name, datatype, array.size(), array_start); - H5Tclose(datatype); + close(datatype); } template <typename DataT, IdType id_type> @@ -215,26 +240,103 @@ class HDF5 DatatypeId datatype = _getDataType<DataT>(); H5Tset_order(datatype, H5T_ORDER_LE); - DatasetId dataset = H5Dopen2(id, name.c_str(), H5P_DEFAULT); + if (not H5Lexists(id, name.c_str(), H5P_DEFAULT)) { + std::ostringstream error_msg; + error_msg << "cannot open dataset '" << name << "'"; + throw NormalError(error_msg.str()); + } - if (dataset < 0) { + DatasetId dataset_id = H5Dopen2(id, name.c_str(), H5P_DEFAULT); + + hid_t dataset_type = H5Dget_type(dataset_id); + if (not H5Tequal(dataset_type, datatype)) { std::ostringstream error_msg; - error_msg << "cannot open dataset \"" << name << "\"\n"; + error_msg << "invalid type for array '" << rang::fgB::yellow << name << rang::fg::reset << "'"; throw NormalError(error_msg.str()); } - const size_t array_size = H5Dget_storage_size(dataset) / sizeof(DataT); - Array<DataT> array{array_size}; + hid_t dataspace_id = H5Dget_space(dataset_id); + int ndims = H5Sget_simple_extent_ndims(dataspace_id); + if (ndims != 1) { + throw NormalError("stored dataset is not an array"); + } + Array<hsize_t> dims(ndims); + H5Sget_simple_extent_dims(dataspace_id, &(dims[0]), NULL); - H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, &(array[0])); + Array<DataT> array{dims[0]}; - H5Dclose(dataset); + H5Dread(dataset_id, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, &(array[0])); - H5Tclose(datatype); + close(dataset_id); + close(datatype); return array; } + template <typename DataT, IdType id_type> + static void + writeTable(const HDFId<id_type>& id, const std::string& name, const Table<DataT>& table) + { + if (H5Lexists(id, name.c_str(), H5P_DEFAULT)) { + std::ostringstream error_msg; + error_msg << "dataset '" << name << "' already exists"; + throw NormalError(error_msg.str()); + } + + static_assert((id_type == IdType::file) or (id_type == IdType::group)); + DatatypeId datatype = _getDataType<DataT>(); + H5Tset_order(datatype, H5T_ORDER_LE); + + const void* table_start = &(table(typename Table<DataT>::index_type(0), typename Table<DataT>::index_type(0))); + + _writeTable(id, name, datatype, table.numberOfRows(), table.numberOfColumns(), table_start); + + close(datatype); + } + + template <typename DataT, IdType id_type> + static Table<DataT> + readTable(const HDFId<id_type>& id, const std::string& name) + { + static_assert((id_type == IdType::file) or (id_type == IdType::group)); + static_assert(not std::is_const_v<DataT>, "DataT must be mutable"); + + DatatypeId datatype = _getDataType<DataT>(); + H5Tset_order(datatype, H5T_ORDER_LE); + + if (not H5Lexists(id, name.c_str(), H5P_DEFAULT)) { + std::ostringstream error_msg; + error_msg << "cannot open dataset '" << name << "'"; + throw NormalError(error_msg.str()); + } + + DatasetId dataset_id = H5Dopen2(id, name.c_str(), H5P_DEFAULT); + + hid_t dataset_type = H5Dget_type(dataset_id); + if (not H5Tequal(dataset_type, datatype)) { + std::ostringstream error_msg; + error_msg << "invalid type for table '" << rang::fgB::yellow << name << rang::fg::reset << "'"; + throw NormalError(error_msg.str()); + } + + hid_t dataspace_id = H5Dget_space(dataset_id); + int ndims = H5Sget_simple_extent_ndims(dataspace_id); + if (ndims != 2) { + throw NormalError("stored dataset is not a table"); + } + Array<hsize_t> dims(ndims); + H5Sget_simple_extent_dims(dataspace_id, &(dims[0]), NULL); + + Table<DataT> table{dims[0], dims[1]}; + + H5Dread(dataset_id, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, &(table(0, 0))); + + close(dataset_id); + close(datatype); + + return table; + } + template <typename DataT, IdType id_type> static void writeAttribute(const HDFId<id_type>& id, const std::string& name, const DataT& value) @@ -246,24 +348,25 @@ class HDF5 const void* v = value.c_str(); _writeAttribute(id, name, datatype, &v); - H5Tclose(datatype); + close(datatype); } else if constexpr (std::is_same_v<std::string_view, std::decay_t<DataT>>) { DatatypeId datatype = H5Tcopy(H5T_C_S1); H5Tset_size(datatype, H5T_VARIABLE); const void* v = &(value[0]); _writeAttribute(id, name, datatype, &v); - H5Tclose(datatype); + close(datatype); } else if constexpr (std::is_same_v<DataT, const char*> or std::is_same_v<DataT, char*>) { DatatypeId datatype = H5Tcopy(H5T_C_S1); H5Tset_size(datatype, H5T_VARIABLE); _writeAttribute(id, name, datatype, &value); - H5Tclose(datatype); - } else if constexpr (std::is_arithmetic_v<std::decay_t<DataT>>) { + close(datatype); + } else if constexpr (std::is_arithmetic_v<std::decay_t<DataT>> or is_tiny_vector_v<std::decay_t<DataT>> or + is_tiny_matrix_v<std::decay_t<DataT>>) { DatatypeId datatype = _getDataType<DataT>(); _writeAttribute(id, name, datatype, &value); - H5Tclose(datatype); + close(datatype); } else { static_assert(std::is_same_v<std::decay_t<DataT>, std::string>, "unexpected attribute type"); } @@ -284,24 +387,49 @@ class HDF5 DataT value = v; H5free_memory(v); - H5Tclose(datatype); + close(datatype); return value; - } else if constexpr (std::is_arithmetic_v<std::decay_t<DataT>>) { + } else if constexpr (std::is_arithmetic_v<std::decay_t<DataT>> or is_tiny_vector_v<DataT> or + is_tiny_matrix_v<DataT>) { DatatypeId datatype = _getDataType<DataT>(); DataT value; _readAttribute(id, name, datatype, &value); - H5Tclose(datatype); + close(datatype); return value; } else { static_assert(std::is_same_v<std::decay_t<DataT>, std::string>, "unexpected attribute type"); } } + static void + close(const AttributeId& attribute_id) + { + H5Aclose(attribute_id); + } + + static void + close(const DatasetId& dataset_id) + { + H5Dclose(dataset_id); + } + + static void + close(const DatatypeId& file_id) + { + H5Tclose(file_id); + } + static void close(const FileId& file_id) { H5Fclose(file_id); } + + static void + close(const GroupId& file_id) + { + H5Gclose(file_id); + } }; #endif // PUGS_HAS_HDF5 diff --git a/tests/test_HDF5.cpp b/tests/test_HDF5.cpp index bb1ceabe5be84cca3223dfc25c53508a8c6749ab..2d424b0635d3af54c9c0978fdef02ecc6af0b1af 100644 --- a/tests/test_HDF5.cpp +++ b/tests/test_HDF5.cpp @@ -3,6 +3,9 @@ #include <utils/HDF5.hpp> +#include <algebra/TinyMatrix.hpp> +#include <algebra/TinyVector.hpp> + #ifdef PUGS_HAS_HDF5 #include <filesystem> @@ -11,19 +14,116 @@ TEST_CASE("HDF5", "[utils]") { + H5Eset_auto(H5E_DEFAULT, NULL, NULL); + static const std::string filename = [&]() -> std::filesystem::path { std::string template_temp_dir = std::filesystem::temp_directory_path() / "test_hdf5_XXXXXX"; return std::string{std::filesystem::path{mkdtemp(&template_temp_dir[0])}} + std::string{".h5"}; }(); + auto is_same_array = [](auto f, auto g) -> bool { + if (f.size() != g.size()) { + return false; + } + for (size_t i = 0; i < f.size(); ++i) { + if (f[i] != g[i]) { + return false; + } + } + return true; + }; + + auto is_same_table = [](auto f, auto g) -> bool { + if ((f.numberOfRows() != g.numberOfRows()) or (f.numberOfColumns() != g.numberOfColumns())) { + return false; + } + for (size_t i = 0; i < f.numberOfRows(); ++i) { + for (size_t j = 0; j < f.numberOfColumns(); ++j) { + if (f(i, j) != g(i, j)) { + return false; + } + } + } + return true; + }; + + Array<const int> int_array = [] { + Array<int> array(12); + for (size_t i = 0; i < array.size(); ++i) { + array[i] = 3 * i + 2; + } + return array; + }(); + + Array<const double> double_array = [] { + Array<double> array(10); + for (size_t i = 0; i < array.size(); ++i) { + array[i] = 3.7 * i - 2.3; + } + return array; + }(); + + Array<const TinyVector<3, double>> R3_array = [] { + Array<TinyVector<3, double>> array(9); + for (size_t i = 0; i < array.size(); ++i) { + array[i] = TinyVector<3, double>{3.7 * i - 2.3, 2.3 + i, -3.2 * i - 1}; + } + return array; + }(); + + Array<const TinyMatrix<2, 3, int>> Z23_array = [] { + Array<TinyMatrix<2, 3, int>> array(13); + for (int i = 0; i < static_cast<int>(array.size()); ++i) { + array[i] = TinyMatrix<2, 3, int>{3 * i - 2, 2 + i, // + -3 * i - 1, 2 + i, // + 2 * i + 3, 3 - 2 * i}; + } + return array; + }(); + + Table<const int> int_table = [] { + Table<int> table(12, 3); + for (size_t i = 0; i < table.numberOfRows(); ++i) { + for (size_t j = 0; j < table.numberOfColumns(); ++j) { + table(i, j) = 3 * i + 2 - 2 * j; + } + } + return table; + }(); + + Table<const double> double_table = [] { + Table<double> table(10, 5); + for (size_t i = 0; i < table.numberOfRows(); ++i) { + for (size_t j = 0; j < table.numberOfColumns(); ++j) { + table(i, j) = 3.2 * i + 2 - 2.5 * j; + } + } + return table; + }(); + SECTION("create file") { + { + std::stringstream error_msg; + error_msg << "error: cannot open file '" << filename << "' for writing"; + REQUIRE_THROWS_WITH(HDF5::openFileRW(filename), error_msg.str()); + } + + { + std::stringstream error_msg; + error_msg << "error: cannot open file '" << filename << "' for reading"; + REQUIRE_THROWS_WITH(HDF5::openFileRO(filename), error_msg.str()); + } + HDF5::FileId file_id = HDF5::create(filename); REQUIRE(file_id >= 0); HDF5::close(file_id); + + REQUIRE_THROWS_WITH(HDF5::create("/a/non/existing/path/to/a/file"), + "error: cannot create file '/a/non/existing/path/to/a/file'"); } - SECTION("open file and add attributes") + SECTION("add attributes at file root") { HDF5::FileId file_id = HDF5::openFileRW(filename); REQUIRE(file_id >= 0); @@ -49,14 +149,123 @@ TEST_CASE("HDF5", "[utils]") HDF5::writeAttribute(file_id, "long double attr", static_cast<long double>(3.14159265)); HDF5::writeAttribute(file_id, "foo attr", std::string{"foobar"}); + HDF5::writeAttribute(file_id, "tinyvector attr", TinyVector<2, int>(2, 3)); + HDF5::writeAttribute(file_id, "tinymatrix attr", TinyMatrix<4, 3, double>(1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12)); + + REQUIRE_THROWS_WITH(HDF5::writeAttribute(file_id, "foo attr", std::string{"bar"}), + "error: attribute 'foo attr' already exists"); + HDF5::close(file_id); } - SECTION("open file ro") + SECTION("add data at file root") { HDF5::FileId file_id = HDF5::openFileRW(filename); REQUIRE(file_id >= 0); + HDF5::writeArray(file_id, "int array", int_array); + HDF5::writeArray(file_id, "double array", double_array); + HDF5::writeArray(file_id, "R3 array", R3_array); + + HDF5::writeTable(file_id, "int table", int_table); + HDF5::writeTable(file_id, "double table", double_table); + + REQUIRE_THROWS_WITH(HDF5::writeArray(file_id, "double array", double_array), + "error: dataset 'double array' already exists"); + + REQUIRE_THROWS_WITH(HDF5::writeTable(file_id, "double table", double_table), + "error: dataset 'double table' already exists"); + + HDF5::close(file_id); + } + + SECTION("add groups and group attributes") + { + HDF5::FileId file_id = HDF5::openFileRW(filename); + REQUIRE(file_id >= 0); + + HDF5::GroupId group_A_id = HDF5::createOrOpenGroup(file_id, "group_A"); + { + REQUIRE_NOTHROW(HDF5::createOrOpenGroup(file_id, "group_A")); + } + REQUIRE(group_A_id >= 0); + + HDF5::writeAttribute(group_A_id, "char attr", static_cast<char>(-'c')); + HDF5::writeAttribute(group_A_id, "uchar attr", static_cast<unsigned char>('a')); + HDF5::writeAttribute(group_A_id, "short attr", static_cast<short>(-23)); + HDF5::writeAttribute(group_A_id, "unsigned short attr", static_cast<unsigned short>(31)); + HDF5::writeAttribute(group_A_id, "int attr", static_cast<int>(-121)); + HDF5::writeAttribute(group_A_id, "unsigned int attr", static_cast<unsigned int>(263)); + HDF5::writeAttribute(group_A_id, "long int attr", static_cast<long int>(-112)); + HDF5::writeAttribute(group_A_id, "long unsigned int attr", static_cast<unsigned long int>(103)); + HDF5::writeAttribute(group_A_id, "int8_t attr", static_cast<int8_t>(-215)); + HDF5::writeAttribute(group_A_id, "uint8_t attr", static_cast<uint8_t>(237)); + HDF5::writeAttribute(group_A_id, "int16_t attr", static_cast<int16_t>(-124)); + HDF5::writeAttribute(group_A_id, "uint16_t attr", static_cast<uint16_t>(244)); + HDF5::writeAttribute(group_A_id, "int32_t attr", static_cast<int32_t>(-1243)); + HDF5::writeAttribute(group_A_id, "uint32_t attr", static_cast<uint32_t>(1427)); + HDF5::writeAttribute(group_A_id, "int64_t attr", static_cast<int64_t>(-15156)); + HDF5::writeAttribute(group_A_id, "uint64_t attr", static_cast<uint64_t>(19999)); + HDF5::writeAttribute(group_A_id, "float attr", static_cast<float>(1.4142)); + HDF5::writeAttribute(group_A_id, "double attr", static_cast<double>(1.4142135)); + HDF5::writeAttribute(group_A_id, "long double attr", static_cast<long double>(1.41421356237)); + HDF5::writeAttribute(group_A_id, "foo attr", std::string{"foobar"}); + + HDF5::writeArray(group_A_id, "int array", int_array); + HDF5::writeTable(group_A_id, "int table", int_table); + + REQUIRE_THROWS_WITH(HDF5::writeArray(group_A_id, "int array", int_array), + "error: dataset 'int array' already exists"); + REQUIRE_THROWS_WITH(HDF5::writeTable(group_A_id, "int table", int_table), + "error: dataset 'int table' already exists"); + + HDF5::GroupId group_B_id = HDF5::createOrOpenGroup(group_A_id, "group_B"); + REQUIRE(group_B_id >= 0); + { + REQUIRE_NOTHROW(HDF5::createOrOpenGroup(group_A_id, "group_B")); + } + + HDF5::writeAttribute(group_B_id, "char attr", static_cast<char>(-'e')); + HDF5::writeAttribute(group_B_id, "uchar attr", static_cast<unsigned char>('s')); + HDF5::writeAttribute(group_B_id, "short attr", static_cast<short>(-133)); + HDF5::writeAttribute(group_B_id, "unsigned short attr", static_cast<unsigned short>(33)); + HDF5::writeAttribute(group_B_id, "int attr", static_cast<int>(-131)); + HDF5::writeAttribute(group_B_id, "unsigned int attr", static_cast<unsigned int>(323)); + HDF5::writeAttribute(group_B_id, "long int attr", static_cast<long int>(-1041)); + HDF5::writeAttribute(group_B_id, "long unsigned int attr", static_cast<unsigned long int>(1301)); + HDF5::writeAttribute(group_B_id, "int8_t attr", static_cast<int8_t>(-325)); + HDF5::writeAttribute(group_B_id, "uint8_t attr", static_cast<uint8_t>(267)); + HDF5::writeAttribute(group_B_id, "int16_t attr", static_cast<int16_t>(-194)); + HDF5::writeAttribute(group_B_id, "uint16_t attr", static_cast<uint16_t>(824)); + HDF5::writeAttribute(group_B_id, "int32_t attr", static_cast<int32_t>(-1234)); + HDF5::writeAttribute(group_B_id, "uint32_t attr", static_cast<uint32_t>(1437)); + HDF5::writeAttribute(group_B_id, "int64_t attr", static_cast<int64_t>(-15415)); + HDF5::writeAttribute(group_B_id, "uint64_t attr", static_cast<uint64_t>(18999)); + HDF5::writeAttribute(group_B_id, "float attr", static_cast<float>(2.718)); + HDF5::writeAttribute(group_B_id, "double attr", static_cast<double>(2.718281)); + HDF5::writeAttribute(group_B_id, "long double attr", static_cast<long double>(2.71828182846)); + HDF5::writeAttribute(group_B_id, "foo attr", std::string{"foobar"}); + + HDF5::writeArray(group_B_id, "double array", double_array); + HDF5::writeArray(group_B_id, "Z23 array", Z23_array); + REQUIRE_THROWS_WITH(HDF5::writeArray(group_B_id, "double array", double_array), + "error: dataset 'double array' already exists"); + + HDF5::close(group_A_id); + HDF5::close(group_B_id); + HDF5::close(file_id); + } + + SECTION("open file ro and open groups") + { + std::filesystem::permissions(filename, + std::filesystem::perms::owner_write | std::filesystem::perms::group_all | + std::filesystem::perms::others_all, + std::filesystem::perm_options::remove); + + HDF5::FileId file_id = HDF5::openFileRO(filename); + REQUIRE(file_id >= 0); + REQUIRE(HDF5::readAttribute<char>(file_id, "char attr") == static_cast<char>(-'b')); REQUIRE(HDF5::readAttribute<unsigned char>(file_id, "uchar attr") == static_cast<unsigned char>('a')); REQUIRE(HDF5::readAttribute<short>(file_id, "short attr") == static_cast<short>(-13)); @@ -79,11 +288,100 @@ TEST_CASE("HDF5", "[utils]") REQUIRE(HDF5::readAttribute<long double>(file_id, "long double attr") == static_cast<long double>(3.14159265)); REQUIRE(HDF5::readAttribute<std::string>(file_id, "foo attr") == std::string{"foobar"}); + REQUIRE(HDF5::readAttribute<TinyVector<2, int>>(file_id, "tinyvector attr") == TinyVector<2, int>(2, 3)); + REQUIRE(HDF5::readAttribute<TinyMatrix<4, 3, double>>(file_id, "tinymatrix attr") == + TinyMatrix<4, 3, double>(1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12)); + + REQUIRE(is_same_array(HDF5::readArray<TinyVector<3, double>>(file_id, "R3 array"), R3_array)); + REQUIRE(is_same_array(HDF5::readArray<int>(file_id, "int array"), int_array)); + REQUIRE(is_same_array(HDF5::readArray<double>(file_id, "double array"), double_array)); + + REQUIRE(is_same_table(HDF5::readTable<int>(file_id, "int table"), int_table)); + REQUIRE(is_same_table(HDF5::readTable<double>(file_id, "double table"), double_table)); + REQUIRE_THROWS_WITH(HDF5::readAttribute<double>(file_id, "invalid attr"), "error: cannot find attribute 'invalid attr'"); REQUIRE_THROWS_WITH(HDF5::readAttribute<double>(file_id, "int64_t attr"), "error: invalid type for attribute 'int64_t attr'"); + REQUIRE_THROWS_WITH(HDF5::readArray<float>(file_id, "invalid array"), "error: cannot open dataset 'invalid array'"); + REQUIRE_THROWS_WITH(HDF5::readArray<float>(file_id, "double array"), + "error: invalid type for array 'double array'"); + + REQUIRE_THROWS_WITH(HDF5::readTable<float>(file_id, "invalid table"), "error: cannot open dataset 'invalid table'"); + REQUIRE_THROWS_WITH(HDF5::readTable<float>(file_id, "double table"), + "error: invalid type for table 'double table'"); + + REQUIRE_THROWS_WITH(HDF5::readArray<double>(file_id, "double table"), "error: stored dataset is not an array"); + REQUIRE_THROWS_WITH(HDF5::readTable<double>(file_id, "double array"), "error: stored dataset is not a table"); + + HDF5::GroupId group_A_id = HDF5::openGroup(file_id, "group_A"); + REQUIRE(group_A_id >= 0); + + REQUIRE(HDF5::readAttribute<char>(group_A_id, "char attr") == static_cast<char>(-'c')); + REQUIRE(HDF5::readAttribute<unsigned char>(group_A_id, "uchar attr") == static_cast<unsigned char>('a')); + REQUIRE(HDF5::readAttribute<short>(group_A_id, "short attr") == static_cast<short>(-23)); + REQUIRE(HDF5::readAttribute<unsigned short>(group_A_id, "unsigned short attr") == static_cast<unsigned short>(31)); + REQUIRE(HDF5::readAttribute<int>(group_A_id, "int attr") == static_cast<int>(-121)); + REQUIRE(HDF5::readAttribute<unsigned int>(group_A_id, "unsigned int attr") == static_cast<unsigned int>(263)); + REQUIRE(HDF5::readAttribute<long int>(group_A_id, "long int attr") == static_cast<long int>(-112)); + REQUIRE(HDF5::readAttribute<unsigned long int>(group_A_id, "long unsigned int attr") == + static_cast<unsigned long int>(103)); + REQUIRE(HDF5::readAttribute<int8_t>(group_A_id, "int8_t attr") == static_cast<int8_t>(-215)); + REQUIRE(HDF5::readAttribute<uint8_t>(group_A_id, "uint8_t attr") == static_cast<uint8_t>(237)); + REQUIRE(HDF5::readAttribute<int16_t>(group_A_id, "int16_t attr") == static_cast<int16_t>(-124)); + REQUIRE(HDF5::readAttribute<uint16_t>(group_A_id, "uint16_t attr") == static_cast<uint16_t>(244)); + REQUIRE(HDF5::readAttribute<int32_t>(group_A_id, "int32_t attr") == static_cast<int32_t>(-1243)); + REQUIRE(HDF5::readAttribute<uint32_t>(group_A_id, "uint32_t attr") == static_cast<uint32_t>(1427)); + REQUIRE(HDF5::readAttribute<int64_t>(group_A_id, "int64_t attr") == static_cast<int64_t>(-15156)); + REQUIRE(HDF5::readAttribute<uint64_t>(group_A_id, "uint64_t attr") == static_cast<uint64_t>(19999)); + REQUIRE(HDF5::readAttribute<float>(group_A_id, "float attr") == static_cast<float>(1.4142)); + REQUIRE(HDF5::readAttribute<double>(group_A_id, "double attr") == static_cast<double>(1.4142135)); + REQUIRE(HDF5::readAttribute<long double>(group_A_id, "long double attr") == + static_cast<long double>(1.41421356237)); + REQUIRE(HDF5::readAttribute<std::string>(group_A_id, "foo attr") == std::string{"foobar"}); + + REQUIRE(is_same_array(HDF5::readArray<int>(group_A_id, "int array"), int_array)); + REQUIRE_THROWS_WITH(HDF5::readArray<double>(group_A_id, "double array"), + "error: cannot open dataset 'double array'"); + + HDF5::GroupId group_B_id = HDF5::openGroup(group_A_id, "group_B"); + REQUIRE(group_B_id >= 0); + + REQUIRE_THROWS_WITH(HDF5::openGroup(file_id, "group_B"), "error: cannot open group 'group_B'"); + REQUIRE_THROWS_WITH(HDF5::openGroup(group_A_id, "group_C"), "error: cannot open group 'group_C'"); + + REQUIRE(HDF5::readAttribute<char>(group_B_id, "char attr") == static_cast<char>(-'e')); + REQUIRE(HDF5::readAttribute<unsigned char>(group_B_id, "uchar attr") == static_cast<unsigned char>('s')); + REQUIRE(HDF5::readAttribute<short>(group_B_id, "short attr") == static_cast<short>(-133)); + REQUIRE(HDF5::readAttribute<unsigned short>(group_B_id, "unsigned short attr") == static_cast<unsigned short>(33)); + REQUIRE(HDF5::readAttribute<int>(group_B_id, "int attr") == static_cast<int>(-131)); + REQUIRE(HDF5::readAttribute<unsigned int>(group_B_id, "unsigned int attr") == static_cast<unsigned int>(323)); + REQUIRE(HDF5::readAttribute<long int>(group_B_id, "long int attr") == static_cast<long int>(-1041)); + REQUIRE(HDF5::readAttribute<unsigned long int>(group_B_id, "long unsigned int attr") == + static_cast<unsigned long int>(1301)); + REQUIRE(HDF5::readAttribute<int8_t>(group_B_id, "int8_t attr") == static_cast<int8_t>(-325)); + REQUIRE(HDF5::readAttribute<uint8_t>(group_B_id, "uint8_t attr") == static_cast<uint8_t>(267)); + REQUIRE(HDF5::readAttribute<int16_t>(group_B_id, "int16_t attr") == static_cast<int16_t>(-194)); + REQUIRE(HDF5::readAttribute<uint16_t>(group_B_id, "uint16_t attr") == static_cast<uint16_t>(824)); + REQUIRE(HDF5::readAttribute<int32_t>(group_B_id, "int32_t attr") == static_cast<int32_t>(-1234)); + REQUIRE(HDF5::readAttribute<uint32_t>(group_B_id, "uint32_t attr") == static_cast<uint32_t>(1437)); + REQUIRE(HDF5::readAttribute<int64_t>(group_B_id, "int64_t attr") == static_cast<int64_t>(-15415)); + REQUIRE(HDF5::readAttribute<uint64_t>(group_B_id, "uint64_t attr") == static_cast<uint64_t>(18999)); + REQUIRE(HDF5::readAttribute<float>(group_B_id, "float attr") == static_cast<float>(2.718)); + REQUIRE(HDF5::readAttribute<double>(group_B_id, "double attr") == static_cast<double>(2.718281)); + REQUIRE(HDF5::readAttribute<long double>(group_B_id, "long double attr") == + static_cast<long double>(2.71828182846)); + + REQUIRE(HDF5::readAttribute<std::string>(group_B_id, "foo attr") == std::string{"foobar"}); + + REQUIRE(is_same_array(HDF5::readArray<double>(group_B_id, "double array"), double_array)); + REQUIRE(is_same_array(HDF5::readArray<TinyMatrix<2, 3, int>>(group_B_id, "Z23 array"), Z23_array)); + + REQUIRE_THROWS_WITH(HDF5::readArray<int>(group_B_id, "int array"), "error: cannot open dataset 'int array'"); + + HDF5::close(group_A_id); + HDF5::close(group_B_id); HDF5::close(file_id); } }