Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Use HDF5 1.10 API #4352

Merged
merged 5 commits into from
Dec 1, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
4 changes: 2 additions & 2 deletions CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -653,7 +653,7 @@ else()
set(HDF5_USE_STATIC_LIBRARIES off)
endif()

find_package(HDF5 COMPONENTS C)
find_package(HDF5 1.10 COMPONENTS C)

if(HDF5_FOUND)
if(HDF5_IS_PARALLEL)
Expand Down Expand Up @@ -694,7 +694,7 @@ if(HDF5_FOUND)

add_library(IO::HDF5 INTERFACE IMPORTED)
target_include_directories(IO::HDF5 INTERFACE "${HDF5_INCLUDE_DIR}")
target_compile_definitions(IO::HDF5 INTERFACE "H5_USE_16_API")
target_compile_definitions(IO::HDF5 INTERFACE "H5_USE_110_API")
target_link_libraries(IO::HDF5 INTERFACE "${HDF5_LIBRARIES}")
if(ENABLE_PHDF5)
target_compile_definitions(IO::HDF5 INTERFACE "ENABLE_PHDF5")
Expand Down
2 changes: 1 addition & 1 deletion README.md
Original file line number Diff line number Diff line change
Expand Up @@ -26,7 +26,7 @@ particular emphasis is placed on code quality and reproducibility.
* CMake v3.17.0 or later, build utility, http://www.cmake.org
* BLAS/LAPACK, numerical library. Use vendor and platform-optimized libraries.
* LibXml2, XML parser, http://xmlsoft.org/
* HDF5, portable I/O library, http://www.hdfgroup.org/HDF5/
* HDF5 v1.10.0 or later, portable I/O library, http://www.hdfgroup.org/HDF5/
* BOOST v1.61.0 or newer, peer-reviewed portable C++ source libraries, http://www.boost.org
* FFTW, FFT library, http://www.fftw.org/
* MPI, parallel library. Optional, but a near requirement for production calculations.
Expand Down
3 changes: 2 additions & 1 deletion src/QMCHamiltonians/SkEstimator.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -128,7 +128,8 @@ void SkEstimator::registerCollectables(std::vector<ObservableHelper>& h5desc, hd
kdims[1] = OHMMS_DIM;
std::string kpath = name_ + "/kpoints";
hid_t k_space = H5Screate_simple(2, kdims, NULL);
hid_t k_set = H5Dcreate(file.getFileID(), kpath.c_str(), H5T_NATIVE_DOUBLE, k_space, H5P_DEFAULT);
hid_t k_set =
H5Dcreate(file.getFileID(), kpath.c_str(), H5T_NATIVE_DOUBLE, k_space, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
hid_t mem_space = H5Screate_simple(2, kdims, NULL);
auto* ptr = &(sourcePtcl->getSimulationCell().getKLists().kpts_cart[0][0]);
herr_t ret = H5Dwrite(k_set, H5T_NATIVE_DOUBLE, mem_space, k_space, H5P_DEFAULT, ptr);
Expand Down
4 changes: 2 additions & 2 deletions src/QMCWaveFunctions/PlaneWave/PWParameterSet.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -76,7 +76,7 @@ bool PWParameterSet::getEigVectorType(hid_t h)
oss << "/" << spinTag << 0;
oss << "/eigenvector";
hsize_t dimTot[4];
hid_t dataset = H5Dopen(h, oss.str().c_str());
hid_t dataset = H5Dopen(h, oss.str().c_str(), H5P_DEFAULT);
hid_t dataspace = H5Dget_space(dataset);
rank = H5Sget_simple_extent_ndims(dataspace);
int status_n = H5Sget_simple_extent_dims(dataspace, dimTot, NULL);
Expand Down Expand Up @@ -184,7 +184,7 @@ void PWParameterSet::checkVersion(hdf_archive& h)
{
if (is_manager())
{
hid_t dataset = H5Dopen(h.getFileID(), "version");
hid_t dataset = H5Dopen(h.getFileID(), "version", H5P_DEFAULT);
hid_t datatype = H5Dget_type(dataset);
H5T_class_t classtype = H5Tget_class(datatype);
H5Tclose(datatype);
Expand Down
2 changes: 1 addition & 1 deletion src/io/hdf/hdf_archive.h
Original file line number Diff line number Diff line change
Expand Up @@ -317,7 +317,7 @@ class hdf_archive
if (Mode[NOIO])
return;
hid_t p = group_id.empty() ? file_id : group_id.top();
herr_t status = H5Gunlink(p, aname.c_str());
herr_t status = H5Ldelete(p, aname.c_str(), H5P_DEFAULT);
}
};

Expand Down
12 changes: 6 additions & 6 deletions src/io/hdf/hdf_stl.h
Original file line number Diff line number Diff line change
Expand Up @@ -96,7 +96,7 @@ struct h5data_proxy<std::string>

inline bool read(data_type& ref, hid_t grp, const std::string& aname, hid_t xfer_plist = H5P_DEFAULT)
{
hid_t dataset = H5Dopen(grp, aname.c_str());
hid_t dataset = H5Dopen(grp, aname.c_str(), H5P_DEFAULT);
if (dataset > -1)
{
hid_t datatype = H5Dget_type(dataset);
Expand Down Expand Up @@ -127,11 +127,11 @@ struct h5data_proxy<std::string>
hsize_t dim = 1;

herr_t ret = -1;
hid_t h1 = H5Dopen(grp, aname.c_str());
hid_t h1 = H5Dopen(grp, aname.c_str(), H5P_DEFAULT);
if (h1 < 0) // missing create one
{
hid_t dataspace = H5Screate_simple(1, &dim, NULL);
hid_t dataset = H5Dcreate(grp, aname.c_str(), str80, dataspace, H5P_DEFAULT);
hid_t dataset = H5Dcreate(grp, aname.c_str(), str80, dataspace, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
ret = H5Dwrite(dataset, str80, H5S_ALL, H5S_ALL, xfer_plist, ref.data());
H5Sclose(dataspace);
H5Dclose(dataset);
Expand All @@ -157,7 +157,7 @@ struct h5data_proxy<std::vector<std::string>>
{
hid_t datatype = H5Tcopy(H5T_C_S1);
H5Tset_size(datatype, H5T_VARIABLE);
hid_t dataset = H5Dopen(grp, aname.c_str());
hid_t dataset = H5Dopen(grp, aname.c_str(), H5P_DEFAULT);
std::vector<char*> char_list;
herr_t ret = -1;
if (dataset > -1)
Expand Down Expand Up @@ -198,12 +198,12 @@ struct h5data_proxy<std::vector<std::string>>
for (int i = 0; i < ref.size(); i++)
char_list.push_back(ref[i].data());

hid_t h1 = H5Dopen(grp, aname.c_str());
hid_t h1 = H5Dopen(grp, aname.c_str(), H5P_DEFAULT);
herr_t ret = -1;
if (h1 < 0) // missing create one
{
hid_t dataspace = H5Screate_simple(1, &dim, NULL);
hid_t dataset = H5Dcreate(grp, aname.c_str(), datatype, dataspace, H5P_DEFAULT);
hid_t dataset = H5Dcreate(grp, aname.c_str(), datatype, dataspace, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
ret = H5Dwrite(dataset, datatype, H5S_ALL, H5S_ALL, xfer_plist, char_list.data());
H5Sclose(dataspace);
H5Dclose(dataset);
Expand Down
26 changes: 13 additions & 13 deletions src/io/hdf/hdf_wrapper_functions.h
Original file line number Diff line number Diff line change
Expand Up @@ -47,7 +47,7 @@ inline bool getDataShape(hid_t grp, const std::string& aname, std::vector<IT>& s
using TSpaceType = h5_space_type<T, 0>;
TSpaceType TSpace;

hid_t h1 = H5Dopen(grp, aname.c_str());
hid_t h1 = H5Dopen(grp, aname.c_str(), H5P_DEFAULT);
hid_t dataspace = H5Dget_space(h1);
int rank = H5Sget_simple_extent_ndims(dataspace);

Expand Down Expand Up @@ -124,7 +124,7 @@ inline bool h5d_read(hid_t grp, const std::string& aname, T* first, hid_t xfer_p
{
if (grp < 0)
return true;
hid_t h1 = H5Dopen(grp, aname.c_str());
hid_t h1 = H5Dopen(grp, aname.c_str(), H5P_DEFAULT);
if (h1 < 0)
return false;
hid_t h5d_type_id = get_h5_datatype(*first);
Expand All @@ -144,12 +144,12 @@ inline bool h5d_write(hid_t grp,
if (grp < 0)
return true;
hid_t h5d_type_id = get_h5_datatype(*first);
hid_t h1 = H5Dopen(grp, aname.c_str());
hid_t h1 = H5Dopen(grp, aname.c_str(), H5P_DEFAULT);
herr_t ret = -1;
if (h1 < 0) //missing create one
{
hid_t dataspace = H5Screate_simple(ndims, dims, NULL);
hid_t dataset = H5Dcreate(grp, aname.c_str(), h5d_type_id, dataspace, H5P_DEFAULT);
hid_t dataset = H5Dcreate(grp, aname.c_str(), h5d_type_id, dataspace, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
ret = H5Dwrite(dataset, h5d_type_id, H5S_ALL, H5S_ALL, xfer_plist, first);
H5Sclose(dataspace);
H5Dclose(dataset);
Expand Down Expand Up @@ -180,7 +180,7 @@ bool h5d_read(hid_t grp,
{
if (grp < 0)
return true;
hid_t h1 = H5Dopen(grp, aname.c_str());
hid_t h1 = H5Dopen(grp, aname.c_str(), H5P_DEFAULT);
if (h1 < 0)
return false;
//herr_t ret = H5Dread(h1, h5d_type_id, H5S_ALL, H5S_ALL, xfer_plist, first);
Expand Down Expand Up @@ -223,15 +223,15 @@ inline bool h5d_write(hid_t grp,
if (grp < 0)
return true;
hid_t h5d_type_id = get_h5_datatype(*first);
hid_t h1 = H5Dopen(grp, aname.c_str());
hid_t h1 = H5Dopen(grp, aname.c_str(), H5P_DEFAULT);
hid_t filespace, memspace;
herr_t ret = -1;

const std::vector<hsize_t> ones(ndims, 1);
if (h1 < 0) //missing create one
{
hid_t dataspace = H5Screate_simple(ndims, gcounts, NULL);
hid_t dataset = H5Dcreate(grp, aname.c_str(), h5d_type_id, dataspace, H5P_DEFAULT);
hid_t dataset = H5Dcreate(grp, aname.c_str(), h5d_type_id, dataspace, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);

hid_t filespace = H5Dget_space(dataset);
ret = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offsets, NULL, ones.data(), counts);
Expand Down Expand Up @@ -276,7 +276,7 @@ bool h5d_read(hid_t grp,
{
if (grp < 0)
return true;
hid_t h1 = H5Dopen(grp, aname.c_str());
hid_t h1 = H5Dopen(grp, aname.c_str(), H5P_DEFAULT);
if (h1 < 0)
return false;

Expand Down Expand Up @@ -320,14 +320,14 @@ inline bool h5d_write(hid_t grp,
<< *(mem_gcounts + 2) << " " << *mem_counts << " " << *(mem_counts + 1) << " " << *(mem_counts + 2) << " "
<< *mem_offsets << " " << *(mem_offsets + 1) << " " << *(mem_offsets + 2) << " " << std::endl;
hid_t h5d_type_id = get_h5_datatype(*first);
hid_t h1 = H5Dopen(grp, aname.c_str());
hid_t h1 = H5Dopen(grp, aname.c_str(), H5P_DEFAULT);
herr_t ret = -1;

const std::vector<hsize_t> ones(std::max(ndims, mem_ndims), 1);
if (h1 < 0) //missing create one
{
hid_t dataspace = H5Screate_simple(ndims, gcounts, NULL);
hid_t dataset = H5Dcreate(grp, aname.c_str(), h5d_type_id, dataspace, H5P_DEFAULT);
hid_t dataset = H5Dcreate(grp, aname.c_str(), h5d_type_id, dataspace, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);

hid_t filespace = H5Dget_space(dataset);
ret = H5Sselect_hyperslab(filespace, H5S_SELECT_SET, offsets, NULL, ones.data(), counts);
Expand Down Expand Up @@ -373,7 +373,7 @@ inline bool h5d_append(hid_t grp,
hid_t h5d_type_id = get_h5_datatype(*first);
hid_t dataspace;
hid_t memspace;
hid_t dataset = H5Dopen(grp, aname.c_str());
hid_t dataset = H5Dopen(grp, aname.c_str(), H5P_DEFAULT);
std::vector<hsize_t> max_dims(ndims);
max_dims[0] = H5S_UNLIMITED;
for (int d = 1; d < ndims; ++d)
Expand All @@ -397,7 +397,7 @@ inline bool h5d_append(hid_t grp,
// set chunk size
hid_t cs = H5Pset_chunk(p, ndims, chunk_dims.data());
// create the dataset
dataset = H5Dcreate2(grp, aname.c_str(), h5d_type_id, dataspace, H5P_DEFAULT, p, H5P_DEFAULT);
dataset = H5Dcreate(grp, aname.c_str(), h5d_type_id, dataspace, H5P_DEFAULT, p, H5P_DEFAULT);
// create memory dataspace, size of current buffer
memspace = H5Screate_simple(ndims, dims, NULL);
// write the data for the first time
Expand Down Expand Up @@ -438,7 +438,7 @@ inline bool h5d_append(hid_t grp,
start[0] = current;
end[0] = start[0] + dims[0];
//extend the dataset (file)
herr_t he = H5Dextend(dataset, end.data());
herr_t he = H5Dset_extent(dataset, end.data());
//get the corresponding dataspace (filespace)
dataspace = H5Dget_space(dataset);
//set the extent
Expand Down