DOLFINx 0.10.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
HDF5Interface.h
1// Copyright (C) 2012 Chris N. Richardson and Garth N. Wells
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#include <array>
10#include <cassert>
11#include <chrono>
12#include <cstdint>
13#include <dolfinx/common/log.h>
14#include <filesystem>
15#include <hdf5.h>
16#include <mpi.h>
17#include <numeric>
18#include <string>
19#include <vector>
20
21namespace dolfinx::io::hdf5
22{
23
25template <typename T>
26hid_t hdf5_type()
27{
28 if constexpr (std::is_same_v<T, float>)
29 return H5T_NATIVE_FLOAT;
30 else if constexpr (std::is_same_v<T, double>)
31 return H5T_NATIVE_DOUBLE;
32 else if constexpr (std::is_same_v<T, std::int32_t>)
33 return H5T_NATIVE_INT32;
34 else if constexpr (std::is_same_v<T, std::uint32_t>)
35 return H5T_NATIVE_UINT32;
36 else if constexpr (std::is_same_v<T, std::int64_t>)
37 return H5T_NATIVE_INT64;
38 else if constexpr (std::is_same_v<T, std::uint64_t>)
39 return H5T_NATIVE_UINT64;
40 else if constexpr (std::is_same_v<T, std::uint8_t>)
41 return H5T_NATIVE_UINT8;
42 else if constexpr (std::is_same_v<T, std::size_t>)
43 {
44 throw std::runtime_error(
45 "Cannot determine size of std::size_t. std::size_t is not the same "
46 "size as long or int.");
47 }
48 else
49 {
50 throw std::runtime_error("Cannot get HDF5 primitive data type. No "
51 "specialised function for this data type.");
52 }
53}
54
60hid_t open_file(MPI_Comm comm, const std::filesystem::path& filename,
61 const std::string& mode, bool use_mpi_io);
62
65void close_file(hid_t handle);
66
69void flush_file(hid_t handle);
70
74std::filesystem::path get_filename(hid_t handle);
75
80bool has_dataset(hid_t handle, const std::string& dataset_path);
81
86hid_t open_dataset(hid_t handle, const std::string& path);
87
92std::vector<std::int64_t> get_dataset_shape(hid_t handle,
93 const std::string& dataset_path);
94
101void set_mpi_atomicity(hid_t handle, bool atomic);
102
107bool get_mpi_atomicity(hid_t handle);
108
112void add_group(hid_t handle, const std::string& dataset_path);
113
124template <typename T>
125void write_dataset(hid_t file_handle, const std::string& dataset_path,
126 const T* data, std::array<std::int64_t, 2> range,
127 const std::vector<int64_t>& global_size, bool use_mpi_io,
128 bool use_chunking)
129{
130 // Data rank
131 const int rank = global_size.size();
132 assert(rank != 0);
133 if (rank > 2)
134 {
135 throw std::runtime_error("Cannot write dataset to HDF5 file"
136 "Only rank 1 and rank 2 dataset are supported");
137 }
138
139 // Get HDF5 data type
140 const hid_t h5type = hdf5::hdf5_type<T>();
141
142 // Hyperslab selection parameters
143 std::vector<hsize_t> count(global_size.begin(), global_size.end());
144 count[0] = range[1] - range[0];
145
146 // Data offsets
147 std::vector<hsize_t> offset(rank, 0);
148 offset[0] = range[0];
149
150 // Dataset dimensions
151 const std::vector<hsize_t> dimsf(global_size.begin(), global_size.end());
152
153 // Create a global data space
154 const hid_t filespace0 = H5Screate_simple(rank, dimsf.data(), nullptr);
155 if (filespace0 == H5I_INVALID_HID)
156 throw std::runtime_error("Failed to create HDF5 data space");
157
158 // Set chunking parameters
159 hid_t chunking_properties;
160 if (use_chunking)
161 {
162 // Set chunk size and limit to 1kB min/1MB max
163 hsize_t chunk_size = dimsf[0] / 2;
164 if (chunk_size > 1048576)
165 chunk_size = 1048576;
166 if (chunk_size < 1024)
167 chunk_size = 1024;
168
169 hsize_t chunk_dims[2] = {chunk_size, dimsf[1]};
170 chunking_properties = H5Pcreate(H5P_DATASET_CREATE);
171 H5Pset_chunk(chunking_properties, rank, chunk_dims);
172 }
173 else
174 chunking_properties = H5P_DEFAULT;
175
176 // Check that group exists and recursively create if required
177 const std::string group_name(dataset_path, 0, dataset_path.rfind('/'));
178 add_group(file_handle, group_name);
179
180 // Create global dataset (using dataset_path)
181 const hid_t dset_id
182 = H5Dcreate2(file_handle, dataset_path.c_str(), h5type, filespace0,
183 H5P_DEFAULT, chunking_properties, H5P_DEFAULT);
184 if (dset_id == H5I_INVALID_HID)
185 throw std::runtime_error("Failed to create HDF5 global dataset.");
186
187 // Close global data space
188 if (herr_t status = H5Sclose(filespace0); status < 0)
189 throw std::runtime_error("Failed to close HDF5 global data space.");
190
191 // Create a local data space
192 const hid_t memspace = H5Screate_simple(rank, count.data(), nullptr);
193 if (memspace == H5I_INVALID_HID)
194 throw std::runtime_error("Failed to create HDF5 local data space.");
195
196 // Create a file dataspace within the global space - a hyperslab
197 const hid_t filespace1 = H5Dget_space(dset_id);
198 herr_t status = H5Sselect_hyperslab(filespace1, H5S_SELECT_SET, offset.data(),
199 nullptr, count.data(), nullptr);
200 if (status < 0)
201 throw std::runtime_error("Failed to create HDF5 dataspace.");
202
203 // Set parallel access
204 const hid_t plist_id = H5Pcreate(H5P_DATASET_XFER);
205 if (use_mpi_io)
206 {
207 if (herr_t status = H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_COLLECTIVE);
208 status < 0)
209 {
210 throw std::runtime_error(
211 "Failed to set HDF5 data transfer property list.");
212 }
213 }
214
215 // Write local dataset into selected hyperslab
216 if (H5Dwrite(dset_id, h5type, memspace, filespace1, plist_id, data) < 0)
217 {
218 throw std::runtime_error(
219 "Failed to write HDF5 local dataset into hyperslab.");
220 }
221
222 if (use_chunking)
223 {
224 // Close chunking properties
225 if (H5Pclose(chunking_properties) < 0)
226 throw std::runtime_error("Failed to close HDF5 chunking properties.");
227 }
228
229 // Close dataset collectively
230 if (H5Dclose(dset_id) < 0)
231 throw std::runtime_error("Failed to close HDF5 dataset.");
232
233 // Close hyperslab
234 if (H5Sclose(filespace1) < 0)
235 throw std::runtime_error("Failed to close HDF5 hyperslab.");
236
237 // Close local dataset
238 if (H5Sclose(memspace) < 0)
239 throw std::runtime_error("Failed to close local HDF5 dataset.");
240
241 // Release file-access template
242 if (H5Pclose(plist_id) < 0)
243 throw std::runtime_error("Failed to release HDF5 file-access template.");
244}
245
255template <typename T>
256std::vector<T> read_dataset(hid_t dset_id, std::array<std::int64_t, 2> range,
257 bool allow_cast)
258{
259 auto timer_start = std::chrono::system_clock::now();
260
261 if (!allow_cast)
262 {
263 // Check that HDF5 dataset type and the type T are the same
264
265 hid_t dtype = H5Dget_type(dset_id);
266 if (dtype == H5I_INVALID_HID)
267 throw std::runtime_error("Failed to get HDF5 data type.");
268 if (htri_t eq = H5Tequal(dtype, hdf5::hdf5_type<T>()); eq < 0)
269 throw std::runtime_error("HDF5 datatype equality test failed.");
270 else if (!eq)
271 {
272 H5Tclose(dtype);
273 throw std::runtime_error("Wrong type for reading from HDF5. Use \"h5ls "
274 "-v\" to inspect the types in the HDF5 file.");
275 }
276 }
277
278 // Open dataspace
279 hid_t dataspace = H5Dget_space(dset_id);
280 if (dataspace == H5I_INVALID_HID)
281 throw std::runtime_error("Failed to open HDF5 data space.");
282
283 // Get rank of data set
284 int rank = H5Sget_simple_extent_ndims(dataspace);
285 if (rank < 1)
286 throw std::runtime_error("Failed to get rank of data space.");
287 else if (rank > 2)
288 spdlog::warn("io::hdf5::read_dataset untested for rank > 2.");
289
290 // Allocate data for shape
291 std::vector<hsize_t> shape(rank);
292
293 // Get size in each dimension
294 if (int ndims = H5Sget_simple_extent_dims(dataspace, shape.data(), nullptr);
295 ndims != rank)
296 {
297 throw std::runtime_error("Failed to get dimensionality of dataspace.");
298 }
299
300 // Hyperslab selection
301 std::vector<hsize_t> offset(rank, 0);
302 std::vector<hsize_t> count = shape;
303 if (range[0] != -1 and range[1] != -1)
304 {
305 offset[0] = range[0];
306 count[0] = range[1] - range[0];
307 }
308 else
309 offset[0] = 0;
310
311 // Select a block in the dataset beginning at offset[], with
312 // size=count[]
313 if (herr_t status
314 = H5Sselect_hyperslab(dataspace, H5S_SELECT_SET, offset.data(), nullptr,
315 count.data(), nullptr);
316 status < 0)
317 {
318 throw std::runtime_error("Failed to select HDF5 hyperslab.");
319 }
320
321 // Create a memory dataspace
322 hid_t memspace = H5Screate_simple(rank, count.data(), nullptr);
323 if (memspace == H5I_INVALID_HID)
324 throw std::runtime_error("Failed to create HDF5 dataspace.");
325
326 // Create local data to read into
327 std::vector<T> data(
328 std::reduce(count.begin(), count.end(), 1, std::multiplies{}));
329
330 // Read data on each process
331 hid_t h5type = hdf5::hdf5_type<T>();
332 if (herr_t status
333 = H5Dread(dset_id, h5type, memspace, dataspace, H5P_DEFAULT, data.data());
334 status < 0)
335 {
336 throw std::runtime_error("Failed to read HDF5 data.");
337 }
338
339 // Close dataspace
340 if (herr_t status = H5Sclose(dataspace); status < 0)
341 throw std::runtime_error("Failed to close HDF5 dataspace.");
342
343 // Close memspace
344 if (herr_t status = H5Sclose(memspace); status < 0)
345 throw std::runtime_error("Failed to close HDF5 memory space.");
346
347 auto timer_end = std::chrono::system_clock::now();
348 std::chrono::duration<double> dt = (timer_end - timer_start);
349 double data_rate = data.size() * sizeof(T) / (1e6 * dt.count());
350 spdlog::info("HDF5 Read data rate: {} MB/s", data_rate);
351
352 return data;
353}
354} // namespace dolfinx::io::hdf5
int rank(MPI_Comm comm)
Return process rank for the communicator.
Definition MPI.cpp:64