13 #include <dolfinx/common/log.h>
35 static hid_t
open_file(MPI_Comm comm,
const std::filesystem::path& filename,
36 const std::string& mode,
const bool use_mpi_io);
62 static void write_dataset(
const hid_t handle,
const std::string& dataset_path,
64 const std::array<std::int64_t, 2>& range,
65 const std::vector<std::int64_t>& global_size,
66 bool use_mpi_io,
bool use_chunking);
78 const std::string& dataset_path,
79 const std::array<std::int64_t, 2>& range);
85 static bool has_dataset(
const hid_t handle,
const std::string& dataset_path);
91 static std::vector<std::int64_t>
112 static void add_group(
const hid_t handle,
const std::string& dataset_path);
115 template <
typename T>
116 static hid_t hdf5_type()
118 throw std::runtime_error(
"Cannot get HDF5 primitive data type. "
119 "No specialised function for this data type");
125 inline hid_t HDF5Interface::hdf5_type<float>()
127 return H5T_NATIVE_FLOAT;
131 inline hid_t HDF5Interface::hdf5_type<double>()
133 return H5T_NATIVE_DOUBLE;
137 inline hid_t HDF5Interface::hdf5_type<int>()
139 return H5T_NATIVE_INT;
143 inline hid_t HDF5Interface::hdf5_type<std::int64_t>()
145 return H5T_NATIVE_INT64;
149 inline hid_t HDF5Interface::hdf5_type<std::size_t>()
151 if (
sizeof(std::size_t) ==
sizeof(
unsigned long))
152 return H5T_NATIVE_ULONG;
153 else if (
sizeof(std::size_t) ==
sizeof(
unsigned int))
154 return H5T_NATIVE_UINT;
156 throw std::runtime_error(
"Cannot determine size of std::size_t. "
157 "std::size_t is not the same size as long or int");
160 template <
typename T>
162 const hid_t file_handle,
const std::string& dataset_path,
const T* data,
163 const std::array<std::int64_t, 2>& range,
164 const std::vector<int64_t>& global_size,
bool use_mpi_io,
bool use_chunking)
167 const int rank = global_size.size();
171 throw std::runtime_error(
"Cannot write dataset to HDF5 file"
172 "Only rank 1 and rank 2 dataset are supported");
176 const hid_t h5type = hdf5_type<T>();
179 std::vector<hsize_t> count(global_size.begin(), global_size.end());
180 count[0] = range[1] - range[0];
183 std::vector<hsize_t> offset(rank, 0);
184 offset[0] = range[0];
187 const std::vector<hsize_t> dimsf(global_size.begin(), global_size.end());
190 const hid_t filespace0 = H5Screate_simple(rank, dimsf.data(),
nullptr);
191 if (filespace0 == HDF5_FAIL)
192 throw std::runtime_error(
"Failed to create HDF5 data space");
195 hid_t chunking_properties;
199 hsize_t chunk_size = dimsf[0] / 2;
200 if (chunk_size > 1048576)
201 chunk_size = 1048576;
202 if (chunk_size < 1024)
205 hsize_t chunk_dims[2] = {chunk_size, dimsf[1]};
206 chunking_properties = H5Pcreate(H5P_DATASET_CREATE);
207 H5Pset_chunk(chunking_properties, rank, chunk_dims);
210 chunking_properties = H5P_DEFAULT;
213 const std::string group_name(dataset_path, 0, dataset_path.rfind(
'/'));
214 add_group(file_handle, group_name);
218 = H5Dcreate2(file_handle, dataset_path.c_str(), h5type, filespace0,
219 H5P_DEFAULT, chunking_properties, H5P_DEFAULT);
220 if (dset_id == HDF5_FAIL)
221 throw std::runtime_error(
"Failed to create HDF5 global dataset.");
227 status = H5Sclose(filespace0);
228 if (status == HDF5_FAIL)
229 throw std::runtime_error(
"Failed to close HDF5 global data space.");
232 const hid_t memspace = H5Screate_simple(rank, count.data(),
nullptr);
233 if (memspace == HDF5_FAIL)
234 throw std::runtime_error(
"Failed to create HDF5 local data space.");
237 const hid_t filespace1 = H5Dget_space(dset_id);
238 status = H5Sselect_hyperslab(filespace1, H5S_SELECT_SET, offset.data(),
239 nullptr, count.data(),
nullptr);
240 if (status == HDF5_FAIL)
241 throw std::runtime_error(
"Failed to create HDF5 dataspace.");
244 const hid_t plist_id = H5Pcreate(H5P_DATASET_XFER);
247 #ifdef H5_HAVE_PARALLEL
248 status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);
249 if (status == HDF5_FAIL)
251 throw std::runtime_error(
252 "Failed to set HDF5 data transfer property list.");
256 throw std::runtime_error(
"HDF5 library has not been configured with MPI");
261 status = H5Dwrite(dset_id, h5type, memspace, filespace1, plist_id, data);
262 if (status == HDF5_FAIL)
264 throw std::runtime_error(
265 "Failed to write HDF5 local dataset into hyperslab.");
271 status = H5Pclose(chunking_properties);
272 if (status == HDF5_FAIL)
273 throw std::runtime_error(
"Failed to close HDF5 chunking properties.");
277 status = H5Dclose(dset_id);
278 if (status == HDF5_FAIL)
279 throw std::runtime_error(
"Failed to close HDF5 dataset.");
282 status = H5Sclose(filespace1);
283 if (status == HDF5_FAIL)
284 throw std::runtime_error(
"Failed to close HDF5 hyperslab.");
287 status = H5Sclose(memspace);
288 if (status == HDF5_FAIL)
289 throw std::runtime_error(
"Failed to close local HDF5 dataset.");
292 status = H5Pclose(plist_id);
293 if (status == HDF5_FAIL)
294 throw std::runtime_error(
"Failed to release HDF5 file-access template.");
297 template <
typename T>
298 inline std::vector<T>
300 const std::string& dataset_path,
301 const std::array<std::int64_t, 2>& range)
303 auto timer_start = std::chrono::system_clock::now();
307 = H5Dopen2(file_handle, dataset_path.c_str(), H5P_DEFAULT);
308 if (dset_id == HDF5_FAIL)
309 throw std::runtime_error(
"Failed to open HDF5 global dataset.");
312 const hid_t dataspace = H5Dget_space(dset_id);
313 if (dataspace == HDF5_FAIL)
314 throw std::runtime_error(
"Failed to open HDF5 data space.");
317 const int rank = H5Sget_simple_extent_ndims(dataspace);
320 LOG(WARNING) <<
"HDF5Interface::read_dataset untested for rank > 2.";
323 std::vector<hsize_t> shape(rank);
326 const int ndims = H5Sget_simple_extent_dims(dataspace, shape.data(),
nullptr);
328 throw std::runtime_error(
"Failed to get dimensionality of dataspace");
331 std::vector<hsize_t> offset(rank, 0);
332 std::vector<hsize_t> count = shape;
333 if (range[0] != -1 and range[1] != -1)
335 offset[0] = range[0];
336 count[0] = range[1] - range[0];
343 [[maybe_unused]] herr_t status = H5Sselect_hyperslab(
344 dataspace, H5S_SELECT_SET, offset.data(),
nullptr, count.data(),
nullptr);
345 if (status == HDF5_FAIL)
346 throw std::runtime_error(
"Failed to select HDF5 hyperslab.");
349 const hid_t memspace = H5Screate_simple(rank, count.data(),
nullptr);
350 if (memspace == HDF5_FAIL)
351 throw std::runtime_error(
"Failed to create HDF5 dataspace.");
355 std::reduce(count.begin(), count.end(), 1, std::multiplies{}));
358 const hid_t h5type = hdf5_type<T>();
360 = H5Dread(dset_id, h5type, memspace, dataspace, H5P_DEFAULT, data.data());
361 if (status == HDF5_FAIL)
362 throw std::runtime_error(
"Failed to read HDF5 data.");
365 status = H5Sclose(dataspace);
366 if (status == HDF5_FAIL)
367 throw std::runtime_error(
"Failed to close HDF5 dataspace.");
370 status = H5Sclose(memspace);
371 if (status == HDF5_FAIL)
372 throw std::runtime_error(
"Failed to close HDF5 memory space.");
375 status = H5Dclose(dset_id);
376 if (status == HDF5_FAIL)
377 throw std::runtime_error(
"Failed to close HDF5 dataset.");
379 auto timer_end = std::chrono::system_clock::now();
380 std::chrono::duration<double> dt = (timer_end - timer_start);
381 double data_rate = data.size() *
sizeof(T) / (1e6 * dt.count());
383 LOG(INFO) <<
"HDF5 Read data rate: " << data_rate <<
"MB/s";
This class provides an interface to some HDF5 functionality.
Definition: HDF5Interface.h:27
static void close_file(const hid_t handle)
Close HDF5 file.
Definition: HDF5Interface.cpp:124
static bool has_dataset(const hid_t handle, const std::string &dataset_path)
Check for existence of dataset in HDF5 file.
Definition: HDF5Interface.cpp:153
static std::filesystem::path get_filename(hid_t handle)
Get filename.
Definition: HDF5Interface.cpp:136
static std::vector< std::int64_t > get_dataset_shape(const hid_t handle, const std::string &dataset_path)
Get dataset shape (size of each dimension)
Definition: HDF5Interface.cpp:203
static hid_t open_file(MPI_Comm comm, const std::filesystem::path &filename, const std::string &mode, const bool use_mpi_io)
Open HDF5 and return file descriptor.
Definition: HDF5Interface.cpp:58
static void write_dataset(const hid_t handle, const std::string &dataset_path, const T *data, const std::array< std::int64_t, 2 > &range, const std::vector< std::int64_t > &global_size, bool use_mpi_io, bool use_chunking)
Write data to existing HDF file as defined by range blocks on each process.
static std::vector< T > read_dataset(const hid_t handle, const std::string &dataset_path, const std::array< std::int64_t, 2 > &range)
Read data from a HDF5 dataset "dataset_path" as defined by range blocks on each process.
static void flush_file(const hid_t handle)
Flush data to file to improve data integrity after interruption.
Definition: HDF5Interface.cpp:130
static bool get_mpi_atomicity(const hid_t handle)
Get MPI atomicity. See https://support.hdfgroup.org/HDF5/doc/RM/RM_H5F.html#File-GetMpiAtomicity and ...
Definition: HDF5Interface.cpp:244
static void set_mpi_atomicity(const hid_t handle, const bool atomic)
Set MPI atomicity. See https://support.hdfgroup.org/HDF5/doc/RM/RM_H5F.html#File-SetMpiAtomicity and ...
Definition: HDF5Interface.cpp:236
int rank(MPI_Comm comm)
Return process rank for the communicator.
Definition: MPI.cpp:76
Support for file IO.
Definition: ADIOS2Writers.h:39