kspaceFirstOrder3D-OMP  1.2
The C++ implementation of the k-wave toolbox for the time-domain simulation of acoustic wave fields in 3D
Hdf5File.h
Go to the documentation of this file.
1 /**
2  * @file Hdf5File.h
3  *
4  * @author Jiri Jaros \n
5  * Faculty of Information Technology\n
6  * Brno University of Technology \n
7  * jarosjir@fit.vutbr.cz
8  *
9  * @brief The header file containing the HDF5 related classes.
10  *
11  * @version kspaceFirstOrder3D 2.16
12  *
13  * @date 27 July 2012, 14:14 (created) \n
14  * 04 September 2017, 10:58 (revised)
15  *
16  *
17  *
18  * @section HDF5File HDF5 File Structure
19  *
20  * The CUDA/C++ code has been designed as a standalone application which is not dependent on MATLAB libraries or a MEX
21  * interface. This is of particular importance when using servers and supercomputers without MATLAB support. For this
22  * reason, simulation data must be transferred between the CUDA/C++ code and MATLAB using external input and output
23  * files. These files are stored using the [Hierarchical Data Format HDF5] (http://www.hdfgroup.org/HDF5/). This is a
24  * data model, library, and file format for storing and managing data. It supports a variety of datatypes, and is
25  * designed for flexible and efficient I/O and for high volume and complex data. The HDF5 technology suite includes
26  * tools and applications for managing, manipulating, viewing, and analysing data in the HDF5 format.
27  *
28  *
29  * Each HDF5 file is a container for storing a variety of scientific data and is composed of two primary types of
30  * objects: groups and datasets. An HDF5 group is a structure containing zero or more HDF5 objects, together with
31  * supporting metadata. A HDF5 group can be seen as a disk folder. An HDF5 dataset is a multidimensional array of data
32  * elements, together with supporting metadata. A HDF5 dataset can be seen as a disk file. Any HDF5 group or dataset
33  * may also have an associated attribute list. An HDF5 attribute is a user-defined HDF5 structure that provides extra
34  * information about an HDF5 object. More information can be obtained from the [HDF5 documentation]
35  * (http://www.hdfgroup.org/HDF5/doc/index.html).
36  *
37  *
38  * kspaceFirstOrder3D-CUDA v1.2 uses the file format introduced in version 1.1. The code is happy to work with both
39  * versions (1.0 and 1.1), however when working with an input file of version 1.0, some features are not supported,
40  * namely the cuboid sensor mask, and <tt>u_non_staggered_raw</tt>. When running from within the actual MATLAB K-Wave
41  * Toolbox, the files will always be generated in version 1.1.
42  *
43  * The HDF5 input file for the CUDA/C++ simulation code contains a file header with brief description of the simulation
44  * stored in string attributes, and the root group <tt>'/'</tt> which stores all the simulation properties in the form
45  * of 3D datasets (a complete list of input datasets is given bellow).
46  * The HDF5 checkpoint file contains the same file header as the input file and the root group <tt>'/'</tt> with a few
47  * datasets keeping the actual simulation state. The HDF5 output file contains a file header with the simulation
48  * description as well as performance statistics, such as the simulation time and memory consumption, stored in string
49  * attributes.
50 
51  * The results of the simulation are stored in the root group <tt>'/'</tt> in the form of 3D datasets. If the linear
52  * sensor mask is used, all output quantities are stored as datasets in the root group. If the cuboid corners sensor
53  * mask is used, the sampled quantities form private groups containing datasets on per cuboid basis.
54  *
55  *\verbatim
56 +----------------------------------------------------------------------------------------------------------------------+
57 | Input File / Checkpoint File Header |
58 +----------------------------------------------------------------------------------------------------------------------+
59 | created_by Short description of the tool that created this file |
60 | creation_date Date when the file was created |
61 | file_description Short description of the content of the file (e.g. simulation name) |
62 | file_type Type of the file (input) |
63 | major_version Major version of the file definition (1) |
64 | minor_version Minor version of the file definition (1) |
65 +----------------------------------------------------------------------------------------------------------------------+
66  \endverbatim
67  *
68  * \verbatim
69 +----------------------------------------------------------------------------------------------------------------------+
70 | Output File Header |
71 +----------------------------------------------------------------------------------------------------------------------+
72 | created_by Short description of the tool that created this file |
73 | creation_date Date when the file was created |
74 | file_description Short description of the content of the file (e.g. simulation name) |
75 | file_type Type of the file (output) |
76 | major_version Major version of the file definition (1) |
77 | minor_version Minor version of the file definition (1) |
78 +----------------------------------------------------------------------------------------------------------------------+
79 | host_names List of hosts (computer names) the simulation was executed on |
80 | number_of_cpu_cores Number of CPU cores used for the simulation |
81 | data_loading_phase_execution_time Time taken to load data from the file |
82 | pre-processing_phase_execution_time Time taken to pre-process data |
83 | simulation_phase_execution_time Time taken to run the simulation |
84 | post-processing_phase_execution_time Time taken to complete the post-processing phase |
85 | total_execution_time Total execution time |
86 | peak_core_memory_in_use Peak memory required per core during the simulation |
87 | total_memory_in_use Total Peak memory in use |
88 +----------------------------------------------------------------------------------------------------------------------+
89 \endverbatim
90  *
91  *
92  *
93  * The input and checkpoint file stores all quantities as three dimensional datasets with dimension sizes designed by
94  * <tt>(Nx, Ny, Nz)</tt>. In order to support scalars and 1D and 2D arrays, the unused dimensions are set to 1. For
95  * example, scalar variables are stored with a dimension size of <tt>(1,1,1)</tt>, 1D vectors oriented in y-direction
96  * are stored with a dimension size of <tt>(1, Ny, 1)</tt>, and so on. If the dataset stores a complex variable, the
97  * real and imaginary parts are stored in an interleaved layout and the lowest used dimension size is doubled (i.e.,
98  * Nx for a 3D matrix, Ny for a 1D vector oriented in the y-direction). The datasets are physically stored in row-major
99  * order (in contrast to column-major order used by MATLAB) using either the <tt>'H5T_IEEE_F32LE'</tt> data type for
100  * floating point datasets or <tt>'H5T_STD_U64LE'</tt> for integer based datasets. All the datasets are store under the
101  * root group.
102  *
103  *
104  * The output file of version 1.0 could only store recorded quantities as 3D datasets under the root group. However,
105  * with version 1.1 and the new cuboid corner sensor mask, the sampled quantities may be laid out as 4D quantities
106  * stored under specific groups. The dimensions are always <tt>(Nx, Ny, Nz, Nt)</tt>, every sampled cuboid is stored as
107  * a distinct dataset and the datasets are grouped under a group named by the quantity stored. This makes the file
108  * clearly readable and easy to parse.
109  *
110  *
111  * In order to enable compression and more efficient data processing, big datasets are not stored as monolithic blocks
112  * but broken into chunks that may be compressed by the ZIP library and stored separately. The chunk size is defined
113  * as follows:
114  *
115  * \li <tt> (1M elements, 1, 1) </tt> in the case of 1D variables - index sensor mask (8MB blocks).
116  * \li <tt> (Nx, Ny, 1) </tt> in the case of 3D variables (one 2D slab).
117  * \li <tt> (Nx, Ny, Nz, 1) </tt> in the case of 4D variables (one time step).
118  * \li <tt> (N_sensor_points, 1, 1) </tt> in the case of the output time series (one time step of
119  * the simulation).
120  *
121  *
122  * All datasets have two attributes that specify the content of the dataset. The <tt>'data_type'</tt> attribute
123  * specifies the data type of the dataset. The admissible values are either <tt>'float'</tt> or <tt>'long'</tt>. The
124  * <tt>'domain_type'</tt> attribute specifies the domain of the dataset. The admissible values are either <tt>'real'
125  * </tt> for the real domain or \c 'complex' for the complex domain. The C++ code reads these attributes and checks
126  * their values.
127  *
128  *
129  *
130  * \verbatim
131 +----------------------------------------------------------------------------------------------------------------------+
132 | Input File Datasets |
133 +----------------------------------------------------------------------------------------------------------------------+
134 | Name Size Data type Domain Type Condition of Presence |
135 +----------------------------------------------------------------------------------------------------------------------+
136 | 1. Simulation Flags |
137 +----------------------------------------------------------------------------------------------------------------------+
138 | ux_source_flag (1, 1, 1) long real |
139 | uy_source_flag (1, 1, 1) long real |
140 | uz_source_flag (1, 1, 1) long real |
141 | p_source_flag (1, 1, 1) long real |
142 | p0_source_flag (1, 1, 1) long real |
143 | transducer_source_flag (1, 1, 1) long real |
144 | nonuniform_grid_flag (1, 1, 1) long real must be set to 0 |
145 | nonlinear_flag (1, 1, 1) long real |
146 | absorbing_flag (1, 1, 1) long real |
147 +----------------------------------------------------------------------------------------------------------------------+
148  2. Grid Properties |
149 +----------------------------------------------------------------------------------------------------------------------+
150 | Nx (1, 1, 1) long real |
151 | Ny (1, 1, 1) long real |
152 | Nz (1, 1, 1) long real |
153 | Nt (1, 1, 1) long real |
154 | dt (1, 1, 1) float real |
155 | dx (1, 1, 1) float real |
156 | dy (1, 1, 1) float real |
157 | dz (1, 1, 1) float real |
158 +----------------------------------------------------------------------------------------------------------------------+
159 | 3. Medium Properties |
160 +----------------------------------------------------------------------------------------------------------------------+
161 | 3.1 Regular Medium Properties |
162 | rho0 (Nx, Ny, Nz) float real heterogenous |
163 | (1, 1, 1) float real homogenous |
164 | rho0_sgx (Nx, Ny, Nz) float real heterogenous |
165 | (1, 1, 1) float real homogenous |
166 | rho0_sgy (Nx, Ny, Nz) float real heterogenous |
167 | (1, 1, 1) float real homogenous |
168 | rho0_sgz (Nx, Ny, Nz) float real heterogenous |
169 | (1, 1, 1) float real homogenous |
170 | c0 (Nx, Ny, Nz) float real heterogenous |
171 | (1, 1, 1) float real homogenous |
172 | c_ref (1, 1, 1) float real |
173 | |
174 | 3.2 Nonlinear Medium Properties (defined if (nonlinear_flag == 1)) |
175 | BonA (Nx, Ny, Nz) float real heterogenous |
176 | (1, 1, 1) float real homogenous |
177 | |
178 | 3.3 Absorbing Medium Properties (defined if (absorbing_flag == 1)) |
179 | alpha_coef (Nx, Ny, Nz) float real heterogenous |
180 | (1, 1, 1) float real homogenous |
181 | alpha_power (1, 1, 1) float real |
182 +----------------------------------------------------------------------------------------------------------------------+
183 | 4. Sensor Variables |
184 +----------------------------------------------------------------------------------------------------------------------+
185 | sensor_mask_type (1, 1, 1) long real file version 1.1 |
186 | (0 = index, 1 = corners) |
187 | sensor_mask_index (Nsens, 1, 1) long real file version 1.0 always, |
188 | file version 1.1 and sensor_mask_type == 0 |
189 | sensor_mask_corners (Ncubes, 6, 1) long real file version 1.1 and sensor_mask_type == 1 |
190 +----------------------------------------------------------------------------------------------------------------------+
191 | 5 Source Properties |
192 +----------------------------------------------------------------------------------------------------------------------+
193 | 5.1 Velocity Source Terms (defined if (ux_source_flag == 1 || uy_source_flag == 1 || uz_source_flag == 1)) |
194 | u_source_mode (1, 1, 1) long real |
195 | u_source_many (1, 1, 1) long real |
196 | u_source_index (Nsrc, 1, 1) long real |
197 | ux_source_input (1, Nt_src, 1) float real u_source_many == 0 |
198 | (Nsrc, Nt_src, 1) float real u_source_many == 1 |
199 | uy_source_input (1, Nt_src, 1) float real u_source_many == 0 |
200 | (Nsrc, Nt_src, 1) float real u_source_many == 1 |
201 | uz_source_input (1, Nt_src, 1) float real u_source_many == 0 |
202 | (Nt_src, Nsrc, 1) float real u_source_many == 1 |
203 | |
204 | 5.2 Pressure Source Terms (defined if (p_source_flag == 1)) |
205 | p_source_mode (1, 1, 1) long real |
206 | p_source_many (1, 1, 1) long real |
207 | p_source_index (Nsrc, 1, 1) long real |
208 | p_source_input (Nsrc, Nt_src, 1) float real p_source_many == 1 |
209 | (1, Nt_src, 1) float real p_source_many == 0 |
210 | |
211 | 5.3 Transducer Source Terms (defined if (transducer_source_flag == 1)) |
212 | u_source_index (Nsrc, 1, 1) long real |
213 | transducer_source_input (Nt_src, 1, 1) float real |
214 | delay_mask (Nsrc, 1, 1) float real |
215 | |
216 | 5.4 IVP Source Terms (defined if ( p0_source_flag == 1)) |
217 | p0_source_input (Nx, Ny, Nz) float real |
218 +----------------------------------------------------------------------------------------------------------------------+
219 | 6. K-space and Shift Variables |
220 +----------------------------------------------------------------------------------------------------------------------+
221 | ddx_k_shift_pos_r (Nx/2 + 1, 1, 1) float complex |
222 | ddx_k_shift_neg_r (Nx/2 + 1, 1, 1) float complex |
223 | ddy_k_shift_pos (1, Ny, 1) float complex |
224 | ddy_k_shift_neg (1, Ny, 1) float complex |
225 | ddz_k_shift_pos (1, 1, Nz) float complex |
226 | ddz_k_shift_neg (1, 1, Nz) float complex |
227 | x_shift_neg_r (Nx/2 + 1, 1, 1) float complex file version 1.1 |
228 | y_shift_neg_r (1, Ny/2 + 1, 1) float complex file version 1.1 |
229 | z_shift_neg_r (1, 1, Nz/2) float complex file version 1.1 |
230 +----------------------------------------------------------------------------------------------------------------------+
231 | 7. PML Variables |
232 +----------------------------------------------------------------------------------------------------------------------+
233 | pml_x_size (1, 1, 1) long real |
234 | pml_y_size (1, 1, 1) long real |
235 | pml_z_size (1, 1, 1) long real |
236 | pml_x_alpha (1, 1, 1) float real |
237 | pml_y_alpha (1, 1, 1) float real |
238 | pml_z_alpha (1, 1, 1) float real |
239 | |
240 | pml_x (Nx, 1, 1) float real |
241 | pml_x_sgx (Nx, 1, 1) float real |
242 | pml_y (1, Ny, 1) float real |
243 | pml_y_sgy (1, Ny, 1) float real |
244 | pml_z (1, 1, Nz) float real |
245 | pml_z_sgz (1, 1, Nz) float real |
246 +----------------------------------------------------------------------------------------------------------------------+
247  \endverbatim
248 
249 \verbatim
250 +----------------------------------------------------------------------------------------------------------------------+
251 | Checkpoint File Datasets |
252 +----------------------------------------------------------------------------------------------------------------------+
253 | Name Size Data type Domain Type Condition of Presence |
254 +----------------------------------------------------------------------------------------------------------------------+
255 | 1. Grid Properties |
256 +----------------------------------------------------------------------------------------------------------------------+
257 | Nx (1, 1, 1) long real |
258 | Ny (1, 1, 1) long real |
259 | Nz (1, 1, 1) long real |
260 | Nt (1, 1, 1) long real |
261 | t_index (1, 1, 1) long real |
262 +----------------------------------------------------------------------------------------------------------------------+
263 | 2. Simulation State |
264 +----------------------------------------------------------------------------------------------------------------------+
265 | p (Nx, Ny, Nz) float real |
266 | ux_sgx (Nx, Ny, Nz) float real |
267 | uy_sgy (Nx, Ny, Nz) float real |
268 | uz_sgz (Nx, Ny, Nz) float real |
269 | rhox (Nx, Ny, Nz) float real |
270 | rhoy (Nx, Ny, Nz) float real |
271 | rhoz (Nx, Ny, Nz) float real |
272 +----------------------------------------------------------------------------------------------------------------------+
273 \endverbatim
274 
275 
276  \verbatim
277 +----------------------------------------------------------------------------------------------------------------------+
278 | Output File Datasets |
279 +----------------------------------------------------------------------------------------------------------------------+
280 | Name Size Data type Domain Type Condition of Presence |
281 +----------------------------------------------------------------------------------------------------------------------+
282 | 1. Simulation Flags |
283 +----------------------------------------------------------------------------------------------------------------------+
284 | ux_source_flag (1, 1, 1) long real |
285 | uy_source_flag (1, 1, 1) long real |
286 | uz_source_flag (1, 1, 1) long real |
287 | p_source_flag (1, 1, 1) long real |
288 | p0_source_flag (1, 1, 1) long real |
289 | transducer_source_flag (1, 1, 1) long real |
290 | nonuniform_grid_flag (1, 1, 1) long real |
291 | nonlinear_flag (1, 1, 1) long real |
292 | absorbing_flag (1, 1, 1) long real |
293 | u_source_mode (1, 1, 1) long real if u_source |
294 | u_source_many (1, 1, 1) long real if u_source |
295 | p_source_mode (1, 1, 1) long real if p_source |
296 | p_source_many (1, 1, 1) long real if p_source |
297 +----------------------------------------------------------------------------------------------------------------------+
298 | 2. Grid Properties |
299 +----------------------------------------------------------------------------------------------------------------------+
300 | Nx (1, 1, 1) long real |
301 | Ny (1, 1, 1) long real |
302 | Nz (1, 1, 1) long real |
303 | Nt (1, 1, 1) long real |
304 | dt (1, 1, 1) float real |
305 | dx (1, 1, 1) float real |
306 | dy (1, 1, 1) float real |
307 | dz (1, 1, 1) float real |
308 +----------------------------------------------------------------------------------------------------------------------+
309 | 3. PML Variables |
310 +----------------------------------------------------------------------------------------------------------------------+
311 | pml_x_size (1, 1, 1) long real |
312 | pml_y_size (1, 1, 1) long real |
313 | pml_z_size (1, 1, 1) long real |
314 | pml_x_alpha (1, 1, 1) float real |
315 | pml_y_alpha (1, 1, 1) float real |
316 | pml_z_alpha (1, 1, 1) float real |
317 | |
318 | pml_x (Nx, 1, 1) float real |
319 | pml_x_sgx (Nx, 1, 1) float real |
320 | pml_y (1, Ny, 1) float real |
321 | pml_y_sgy (1, Ny, 1) float real |
322 | pml_z (1, 1, Nz) float real |
323 | pml_z_sgz (1, 1, Nz) float real |
324 +----------------------------------------------------------------------------------------------------------------------+
325 | 4. Sensor Variables (present if --copy_sensor_mask) |
326 +----------------------------------------------------------------------------------------------------------------------+
327 | sensor_mask_type (1, 1, 1) long real file version 1.1 and --copy_sensor_mask |
328 | sensor_mask_index (Nsens, 1, 1) long real file version 1.1 and sensor_mask_type == 0 |
329 | sensor_mask_corners (Ncubes, 6, 1) long real file version 1.1 and sensor_mask_type == 1 |
330 +----------------------------------------------------------------------------------------------------------------------+
331 | 5a. Simulation Results: if sensor_mask_type == 0 (index), or File version == 1.0 |
332 +----------------------------------------------------------------------------------------------------------------------+
333 | p (Nsens, Nt - s, 1) float real -p or --p_raw |
334 | p_rms (Nsens, 1, 1) float real --p_rms |
335 | p_max (Nsens, 1, 1) float real --p_max |
336 | p_min (Nsens, 1, 1) float real --p_min |
337 | p_max_all (Nx, Ny, Nz) float real --p_max_all |
338 | p_min_all (Nx, Ny, Nz) float real --p_min_all |
339 | p_final (Nx, Ny, Nz) float real --p_final |
340 | |
341 | ux (Nsens, Nt - s, 1) float real -u or --u_raw |
342 | uy (Nsens, Nt - s, 1) float real -u or --u_raw |
343 | uz (Nsens, Nt - s, 1) float real -u or --u_raw |
344 | |
345 | ux_non_staggered (Nsens, Nt - s, 1) float real --u_non_staggered_raw (File version ==1.1) |
346 | uy_non_staggered (Nsens, Nt - s, 1) float real --u_non_staggered_raw (File version ==1.1) |
347 | uz_non_staggered (Nsens, Nt - s, 1) float real --u_non_staggered_raw (File version ==1.1) |
348 | |
349 | ux_rms (Nsens, 1, 1) float real --u_rms |
350 | uy_rms (Nsens, 1, 1) float real --u_rms |
351 | uz_rms (Nsens, 1, 1) float real --u_rms |
352 | |
353 | ux_max (Nsens, 1, 1) float real --u_max |
354 | uy_max (Nsens, 1, 1) float real --u_max |
355 | uz_max (Nsens, 1, 1) float real --u_max |
356 | |
357 | ux_min (Nsens, 1, 1) float real --u_min |
358 | uy_min (Nsens, 1, 1) float real --u_min |
359 | uz_min (Nsens, 1, 1) float real --u_min |
360 | |
361 | ux_max_all (Nx, Ny, Nz) float real --u_max_all |
362 | uy_max_all (Nx, Ny, Nz) float real --u_max_all |
363 | uz_max_all (Nx, Ny, Nz) float real --u_max_all |
364 | |
365 | ux_min_all (Nx, Ny, Nz) float real --u_min_all |
366 | uy_min_all (Nx, Ny, Nz) float real --u_min_all |
367 | uz_min_all (Nx, Ny, Nz) float real --u_min_all |
368 | |
369 | ux_final (Nx, Ny, Nz) float real --u_final |
370 | uy_final (Nx, Ny, Nz) float real --u_final |
371 | uz_final (Nx, Ny, Nz) float real --u_final |
372 +----------------------------------------------------------------------------------------------------------------------+
373 | 5b. Simulation Results: if sensor_mask_type == 1 (corners) and file version == 1.1 |
374 +----------------------------------------------------------------------------------------------------------------------+
375 | /p group of datasets, one per cuboid -p or --p_raw |
376 | /p/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
377 | /p/2 (Cx, Cy, Cz, Nt-s) float real 2nd sampled cuboid, etc. |
378 | |
379 | /p_rms group of datasets, one per cuboid --p_rms |
380 | /p_rms/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
381 | |
382 | /p_max group of datasets, one per cuboid --p_max |
383 | /p_max/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
384 | |
385 | /p_min group of datasets, one per cuboid --p_min |
386 | /p_min/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
387 | |
388 | p_max_all (Nx, Ny, Nz) float real --p_max_all |
389 | p_min_all (Nx, Ny, Nz) float real --p_min_all |
390 | p_final (Nx, Ny, Nz) float real --p_final |
391 | |
392 | |
393 | /ux group of datasets, one per cuboid -u or --u_raw |
394 | /ux/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
395 | /uy group of datasets, one per cuboid -u or --u_raw |
396 | /uy/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
397 | /uz group of datasets, one per cuboid -u or --u_raw |
398 | /uz/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
399 | |
400 | /ux_non_staggered group of datasets, one per cuboid --u_non_staggered_raw |
401 | /ux_non_staggered/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
402 | /uy_non_staggered group of datasets, one per cuboid --u_non_staggered_raw |
403 | /uy_non_staggered/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
404 | /uz_non_staggered group of datasets, one per cuboid --u_non_staggered_raw |
405 | /uz_non_staggered/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
406 | |
407 | /ux_rms group of datasets, one per cuboid --u_rms |
408 | /ux_rms/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
409 | /uy_rms group of datasets, one per cuboid --u_rms |
410 | /uy_rms/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
411 | /uz_rms group of datasets, one per cuboid --u_rms |
412 | /uy_rms/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
413 | |
414 | /ux_max group of datasets, one per cuboid --u_max |
415 | /ux_max/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
416 | /uy_max group of datasets, one per cuboid --u_max |
417 | /ux_max/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
418 | /uz_max group of datasets, one per cuboid --u_max |
419 | /ux_max/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
420 | |
421 | /ux_min group of datasets, one per cuboid --u_min |
422 | /ux_min/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
423 | /uy_min group of datasets, one per cuboid --u_min |
424 | /ux_min/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
425 | /uz_min group of datasets, one per cuboid --u_min |
426 | /ux_min/1 (Cx, Cy, Cz, Nt-s) float real 1st sampled cuboid |
427 | |
428 | ux_max_all (Nx, Ny, Nz) float real --u_max_all |
429 | uy_max_all (Nx, Ny, Nz) float real --u_max_all |
430 | uz_max_all (Nx, Ny, Nz) float real --u_max_all |
431 | |
432 | ux_min_all (Nx, Ny, Nz) float real --u_min_all |
433 | uy_min_all (Nx, Ny, Nz) float real --u_min_all |
434 | uz_min_all (Nx, Ny, Nz) float real --u_min_all |
435 | |
436 | ux_final (Nx, Ny, Nz) float real --u_final |
437 | uy_final (Nx, Ny, Nz) float real --u_final |
438 | uz_final (Nx, Ny, Nz) float real --u_final |
439 +----------------------------------------------------------------------------------------------------------------------+
440 \endverbatim
441  *
442  *
443  *
444  * @copyright Copyright (C) 2017 Jiri Jaros and Bradley Treeby.
445  *
446  * This file is part of the C++ extension of the [k-Wave Toolbox](http://www.k-wave.org).
447  *
448  * This file is part of the k-Wave. k-Wave is free software: you can redistribute it and/or modify it under the terms
449  * of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the
450  * License, or (at your option) any later version.
451  *
452  * k-Wave is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied
453  * warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
454  * more details.
455  *
456  * You should have received a copy of the GNU Lesser General Public License along with k-Wave.
457  * If not, see [http://www.gnu.org/licenses/](http://www.gnu.org/licenses/).
458  */
459 
460 
461 #ifndef THDF5_FILE_H
462 #define THDF5_FILE_H
463 
464 
465 #include <hdf5.h>
466 #include <hdf5_hl.h>
467 #include <string>
468 #include <map>
469 
470 // Linux build
471 #ifdef __linux__
472  #include <unistd.h>
473 #endif
474 
475 #ifdef _WIN64
476  #include <io.h>
477 #endif
478 
479 #include <Utils/MatrixNames.h>
480 #include <Utils/DimensionSizes.h>
481 
482 
483 /**
484  * @class Hdf5File
485  * @brief Class wrapping the HDF5 routines.
486  *
487  * This class is responsible for working with HDF5 files. It offers routines to manage files (create, open, close)
488  * as well as creating, reading and modifying the contents (groups and datasets).
489  */
490 class Hdf5File
491 {
492  public:
493 
494  /**
495  * @enum MatrixDataType
496  * @brief HDF5 matrix data type (float or uint64).
497  * @details HDF5 matrix data type (float or uint64).
498  */
499  enum class MatrixDataType
500  {
501  /// The matrix is stored in floating point 32b wide format.
502  kFloat = 0,
503  /// The matrix is stored in fixed point point 64b wide format.
504  kLong = 1
505  };
506 
507  /**
508  * @enum MatrixDomainType
509  * @brief HDF5 Matrix domain type (real or complex).
510  * @details HDF5 Matrix domain type (real or complex).
511  */
512  enum class MatrixDomainType
513  {
514  /// The matrix is defined on real domain.
515  kReal = 0,
516  /// The matrix is defined on complex domain.
517  kComplex = 1
518  };
519 
520  /// Constructor of the class.
521  Hdf5File();
522  /// Copy constructor is not allowed.
523  Hdf5File(const Hdf5File&) = delete;
524  /// Destructor.
525  virtual ~Hdf5File();
526 
527  /// Operator = is not allowed.
528  Hdf5File& operator=(const Hdf5File&) = delete;
529 
530  //------------------------------------------- Basic file operations ----------------------------------------------//
531  /**
532  * @brief Create the HDF5 file.
533  *
534  * The file is always created in read write mode mode, by default overwriting an old file of the same name.
535  * Other HDF5 flags are set to default.
536  *
537  * @param [in] fileName - File name.
538  * @param [in] flags - How to create the file, by default overwrite existing file.
539  * @throw ios:failure - If error happened (file is open or cannot be created).
540  */
541  void create(const std::string& fileName,
542  unsigned int flags = H5F_ACC_TRUNC);
543  /**
544  * @brief Open the HDF5 file.
545  *
546  * The file is opened in read only mode by default. Other HDF5 flags are set to default.
547  *
548  * @param [in] fileName - File name
549  * @param [in] flags - Open mode, by default read only.
550  * @throw ios:failure - If error happened (file not found, file is not a HDF5 file, file is already open).
551  *
552  */
553  void open(const std::string& fileName,
554  unsigned int flags = H5F_ACC_RDONLY);
555  /**
556  * @brief Is the file opened?
557  * @details Is the file opened?
558  * @return true - If the file is opened.
559  */
560  bool isOpen() const { return mFile != H5I_BADID; };
561  /**
562  * @brief Can I access the file.
563  * @details Can the code access the file, e.g. does it exist, do we have enough privileges, etc.
564  *
565  * @param [in] fileName - Name of the file.
566  * @return true - If it is possible to access the file.
567  */
568  static bool canAccess(const std::string& fileName);
569  /**
570  * @brief Close the HDF5 file.
571  * @throw ios::failure - If an error happens
572  */
573  void close();
574 
575  //--------------------------------------------- Group manipulators -----------------------------------------------//
576  /**
577  * @brief Create a HDF5 group at a specified place in the file tree.
578  * @details Other HDF5 flags are set to default.
579  *
580  * @param [in] parentGroup - Where to link the group at.
581  * @param [in] groupName - Group name.
582  * @return A handle to the new group.
583  * @throw ios::failure - If error happens.
584  */
585  hid_t createGroup(const hid_t parentGroup,
586  MatrixName& groupName);
587  /**
588  * @brief Open a HDF5 group at a specified place in the file tree.
589  * @details Other HDF5 flags are set to default.
590  *
591  * @param [in] parentGroup - Parent group.
592  * @param [in] groupName - Group name.
593  * @return A handle to the group.
594  * @throw ios::failure - If error happens.
595  */
596  hid_t openGroup(const hid_t parentGroup,
597  MatrixName& groupName);
598  /**
599  * @brief Close a group.
600  * @param[in] group - Group to close.
601  */
602  void closeGroup(const hid_t group);
603  /**
604  * @brief Get handle to the root group of the file.
605  * @details Get handle to the root group of the file.
606  * @return Handle to the root group.
607  */
608  hid_t getRootGroup() const { return mFile; };
609 
610 
611  //-------------------------------------------- Dataset manipulators ----------------------------------------------//
612  /**
613  * @brief Open a dataset at a specified place in the file tree.
614  * @details Other HDF5 flags are set to default.
615  *
616  * @param [in] parentGroup - Parent group id (can be the file id for root)..
617  * @param [in] datasetName - Dataset name.
618  * @return A handle to open dataset.
619  * @throw ios::failure - If error happens.
620  */
621  hid_t openDataset(const hid_t parentGroup,
622  MatrixName& datasetName);
623 
624  /**
625  * @brief Create a float HDF5 dataset at a specified place in the file tree (3D/4D).
626  * @details Other HDF5 flags are set to default.
627  *
628  * @param [in] parentGroup - Parent group id.
629  * @param [in] datasetName - Dataset name.
630  * @param [in] dimensionSizes - Dimension sizes.
631  * @param [in] chunkSizes - Chunk sizes.
632  * @param [in] matrixDataType - Matrix data type
633  * @param [in] compressionLevel - Compression level.
634  * @return A handle to the new dataset.
635  * @throw ios::failure - If error happens.
636  */
637  hid_t createDataset(const hid_t parentGroup,
638  MatrixName& datasetName,
639  const DimensionSizes& dimensionSizes,
640  const DimensionSizes& chunkSizes,
641  const MatrixDataType matrixDataType,
642  const size_t compressionLevel);
643 
644  /**
645  * @brief Close dataset.
646  * @param [in] dataset - Dataset to close.
647  */
648  void closeDataset(const hid_t dataset);
649 
650 
651  //---------------------------------------- Dataset Read/Write operations -----------------------------------------//
652  /**
653  * @brief Write a hyperslab into the dataset.
654  *
655  * @tparam T - Data type to be written.
656  * @param [in] dataset - Dataset id.
657  * @param [in] position - Position in the dataset.
658  * @param [in] size - Size of the hyperslab.
659  * @param [in] data - Data to be written.
660  * @throw ios::failure - If error happens.
661  * @warning Limited to float and size_t data types.
662  */
663  template<class T>
664  void writeHyperSlab(const hid_t dataset,
665  const DimensionSizes& position,
666  const DimensionSizes& size,
667  const T* data);
668 
669  /**
670  * @brief Write a cuboid selected within the matrixData into a hyperslab.
671  * @details The routine writes 3D cuboid into a 4D dataset (only intended for raw time series).
672  *
673  * @param [in] dataset - Dataset to write MatrixData into.
674  * @param [in] hyperslabPosition - Position in the dataset (hyperslab) - may be 3D/4D.
675  * @param [in] cuboidPosition - Position of the cuboid in MatrixData (what to sample) - must be 3D.
676  * @param [in] cuboidSize - Cuboid size (size of data being sampled) - must by 3D.
677  * @param [in] matrixDimensions - Size of the original matrix (the sampled one).
678  * @param [in] matrixData - C array of matrix data.
679  * @throw ios::failure - If error happens.
680  */
681  void writeCuboidToHyperSlab(const hid_t dataset,
682  const DimensionSizes& hyperslabPosition,
683  const DimensionSizes& cuboidPosition,
684  const DimensionSizes& cuboidSize,
685  const DimensionSizes& matrixDimensions,
686  const float* matrixData);
687 
688  /**
689  * @brief Write sensor data selected by the sensor mask - Occasionally very slow, do not use!
690  *
691  * Write sensor data selected by the sensor mask. The routine picks elements from the MatixData based on
692  * the Sensor Data and store them into a single hyperslab of size [Nsens, 1, 1]
693  *
694  * @param [in] dataset - Dataset to write MaatrixData into
695  * @param [in] hyperslabPosition - 3D position in the dataset (hyperslab)
696  * @param [in] indexSensorSize - Size of the index based sensor mask
697  * @param [in] indexSensorData - Index based sensor mask
698  * @param [in] matrixDimensions - Size of the sampled matrix
699  * @param [in] matrixData - Matrix data
700  * @warning - very slow at this version of HDF5 for orthogonal planes-> DO NOT USE
701  */
702  void writeSensorByMaskToHyperSlab(const hid_t dataset,
703  const DimensionSizes& hyperslabPosition,
704  const size_t indexSensorSize,
705  const size_t* indexSensorData,
706  const DimensionSizes& matrixDimensions,
707  const float* matrixData);
708 
709  /**
710  * @brief Write the scalar value under a specified group
711  * @details No chunks and no compression is used.
712  *
713  * @tparam T - Data type to be written.
714  * @param [in] parentGroup - Where to link the scalar dataset.
715  * @param [in] datasetName - HDF5 dataset name.
716  * @param [in] value - data to be written.
717  * @throw ios::failure - If error happens.
718  * @warning Limited to float and size_t data types
719  */
720  template<class T>
721  void writeScalarValue(const hid_t parentGroup,
722  MatrixName& datasetName,
723  const T value);
724  /**
725  * @brief Read the scalar value under a specified group.
726  *
727  * @tparam T - Data type to be written.
728  * @param [in] parentGroup - Where to link the scalar dataset.
729  * @param [in] datasetName - HDF5 dataset name.
730  * @param [out] value - Data to be read.
731  * @throw ios::failure - If error happens.
732  * @warning Limited to float and size_t data types
733  */
734  template<class T>
735  void readScalarValue(const hid_t parentGroup,
736  MatrixName& datasetName,
737  T& value);
738 
739  /**
740  * @brief Read data from the dataset at a specified place in the file tree.
741  *
742  * @tparam T - Data type to be written.
743  * @param [in] parentGroup - Where is the dataset situated.
744  * @param [in] datasetName - Dataset name.
745  * @param [in] dimensionSizes - Dimension sizes.
746  * @param [out] data - Pointer to data.
747  * @throw ios::failure - If error happens.
748  */
749  template<class T>
750  void readCompleteDataset(const hid_t parentGroup,
751  MatrixName& datasetName,
752  const DimensionSizes& dimensionSizes,
753  T* data);
754 
755 
756  //--------------------------------------- Attributes Read/Write operations ---------------------------------------//
757  /**
758  * @brief Get dimension sizes of the dataset under a specified group.
759  *
760  * @param [in] parentGroup - Where the dataset is.
761  * @param [in] datasetName - Dataset name.
762  * @return Dimension sizes of the dataset.
763  * @throw ios::failure - If error happens.
764  */
765  DimensionSizes getDatasetDimensionSizes(const hid_t parentGroup,
766  MatrixName& datasetName);
767  /**
768  * @brief Get number of dimensions of the dataset under a specified group.
769  *
770  * @param [in] parentGroup - Where the dataset is.
771  * @param [in] datasetName - Dataset name.
772  * @return Number of dimensions.
773  * @throw ios::failure - If error happens.
774  */
775  size_t getDatasetNumberOfDimensions(const hid_t parentGroup,
776  MatrixName& datasetName);
777  /**
778  * @brief Get dataset element count at a specified place in the file tree.
779  *
780  * @param [in] parentGroup - Where the dataset is.
781  * @param [in] datasetName - Dataset name.
782  * @return Number of elements.
783  * @throw ios::failure - If error happens.
784  */
785  size_t getDatasetSize(const hid_t parentGroup,
786  MatrixName& datasetName);
787  /**
788  * @brief Write matrix data type into the dataset at a specified place in the file tree.
789  *
790  * @param [in] parentGroup - Where the dataset is.
791  * @param [in] datasetName - Dataset name.
792  * @param [in] matrixDataType - Matrix data type in the file.
793  * @throw ios::failure - If error happens.
794  */
795  void writeMatrixDataType (const hid_t parentGroup,
796  MatrixName& datasetName,
797  const MatrixDataType& matrixDataType);
798  /**
799  * @brief Write matrix data type into the dataset at a specified place in the file tree.
800  *
801  * @param [in] parentGroup - Where the dataset is.
802  * @param [in] datasetName - Dataset name.
803  * @param [in] matrixDomainType - Matrix domain type.
804  * @throw ios::failure - If error happens.
805  */
806  void writeMatrixDomainType(const hid_t parentGroup,
807  MatrixName& datasetName,
808  const MatrixDomainType& matrixDomainType);
809  /**
810  * @brief Read matrix data type from the dataset at a specified place in the file tree.
811  *
812  * @param [in] parentGroup - Where the dataset is.
813  * @param [in] datasetName - Dataset name.
814  * @return Matrix data type.
815  * @throw ios::failure - If error happens.
816  */
817  MatrixDataType readMatrixDataType(const hid_t parentGroup,
818  MatrixName& datasetName);
819  /**
820  * @brief Read matrix dataset domain type at a specified place in the file tree.
821  *
822  * @param [in] parentGroup - Where the dataset is.
823  * @param [in] datasetName - Dataset name.
824  * @return Matrix domain type.
825  * @throw ios::failure - If error happens.
826  */
827  MatrixDomainType readMatrixDomainType(const hid_t parentGroup,
828  MatrixName& datasetName);
829  /**
830  * @brief Write string attribute into the dataset under the root group.
831  *
832  * @param [in] parentGroup - Where the dataset is.
833  * @param [in] datasetName - Dataset name.
834  * @param [in] attributeName - Attribute name.
835  * @param [in] value - Data to write.
836  * @throw ios::failure - If error happens.
837  */
838  void writeStringAttribute(const hid_t parentGroup,
839  MatrixName& datasetName,
840  MatrixName& attributeName,
841  const std::string& value);
842  /**
843  * @brief Read string attribute from the dataset under the root group.
844  *
845  * @param [in] parentGroup - Where the dataset is.
846  * @param [in] datasetName - Dataset name.
847  * @param [in] attributeName - Attribute name.
848  * @return Attribute value.
849  * @throw ios::failure - If error happens.
850  */
851  std::string readStringAttribute(const hid_t parentGroup,
852  MatrixName& datasetName,
853  MatrixName& attributeName);
854 
855  private:
856 
857  /// String representation of the Domain type in the HDF5 file.
858  static const std::string kMatrixDomainTypeName;
859  /// String representation of the Data type in the HDF5 file.
860  static const std::string kMatrixDataTypeName;
861 
862  /// String representation of different domain types.
863  static const std::string kMatrixDomainTypeNames[];
864  /// String representation of different data types.
865  static const std::string kMatrixDataTypeNames[];
866 
867  /// HDF file handle.
868  hid_t mFile;
869  /// File name.
870  std::string mFileName;
871 }; // Hdf5File
872 //----------------------------------------------------------------------------------------------------------------------
873 
874 
875 #endif /* THDF5_FILE_H */
void writeStringAttribute(const hid_t parentGroup, MatrixName &datasetName, MatrixName &attributeName, const std::string &value)
Write string attribute into the dataset under the root group.
Definition: Hdf5File.cpp:876
hid_t openDataset(const hid_t parentGroup, MatrixName &datasetName)
Open a dataset at a specified place in the file tree.
Definition: Hdf5File.cpp:223
void create(const std::string &fileName, unsigned int flags=H5F_ACC_TRUNC)
Create the HDF5 file.
Definition: Hdf5File.cpp:93
void writeScalarValue(const hid_t parentGroup, MatrixName &datasetName, const T value)
Write the scalar value under a specified group.
Definition: Hdf5File.cpp:559
static const std::string kMatrixDomainTypeName
String representation of the Domain type in the HDF5 file.
Definition: Hdf5File.h:858
hid_t getRootGroup() const
Get handle to the root group of the file.
Definition: Hdf5File.h:608
virtual ~Hdf5File()
Destructor.
Definition: Hdf5File.cpp:83
hid_t openGroup(const hid_t parentGroup, MatrixName &groupName)
Open a HDF5 group at a specified place in the file tree.
Definition: Hdf5File.cpp:196
void open(const std::string &fileName, unsigned int flags=H5F_ACC_RDONLY)
Open the HDF5 file.
Definition: Hdf5File.cpp:117
void writeCuboidToHyperSlab(const hid_t dataset, const DimensionSizes &hyperslabPosition, const DimensionSizes &cuboidPosition, const DimensionSizes &cuboidSize, const DimensionSizes &matrixDimensions, const float *matrixData)
Write a cuboid selected within the matrixData into a hyperslab.
Definition: Hdf5File.cpp:431
size_t getDatasetNumberOfDimensions(const hid_t parentGroup, MatrixName &datasetName)
Get number of dimensions of the dataset under a specified group.
Definition: Hdf5File.cpp:757
DimensionSizes getDatasetDimensionSizes(const hid_t parentGroup, MatrixName &datasetName)
Get dimension sizes of the dataset under a specified group.
Definition: Hdf5File.cpp:731
std::string mFileName
File name.
Definition: Hdf5File.h:870
hid_t mFile
HDF file handle.
Definition: Hdf5File.h:868
static bool canAccess(const std::string &fileName)
Can I access the file.
Definition: Hdf5File.cpp:146
void writeSensorByMaskToHyperSlab(const hid_t dataset, const DimensionSizes &hyperslabPosition, const size_t indexSensorSize, const size_t *indexSensorData, const DimensionSizes &matrixDimensions, const float *matrixData)
Write sensor data selected by the sensor mask - Occasionally very slow, do not use! ...
Definition: Hdf5File.cpp:492
void writeMatrixDomainType(const hid_t parentGroup, MatrixName &datasetName, const MatrixDomainType &matrixDomainType)
Write matrix data type into the dataset at a specified place in the file tree.
Definition: Hdf5File.cpp:807
void closeDataset(const hid_t dataset)
Close dataset.
Definition: Hdf5File.cpp:325
void readScalarValue(const hid_t parentGroup, MatrixName &datasetName, T &value)
Read the scalar value under a specified group.
Definition: Hdf5File.cpp:646
size_t getDatasetSize(const hid_t parentGroup, MatrixName &datasetName)
Get dataset element count at a specified place in the file tree.
Definition: Hdf5File.cpp:775
void writeMatrixDataType(const hid_t parentGroup, MatrixName &datasetName, const MatrixDataType &matrixDataType)
Write matrix data type into the dataset at a specified place in the file tree.
Definition: Hdf5File.cpp:793
std::string readStringAttribute(const hid_t parentGroup, MatrixName &datasetName, MatrixName &attributeName)
Read string attribute from the dataset under the root group.
Definition: Hdf5File.cpp:892
hid_t createGroup(const hid_t parentGroup, MatrixName &groupName)
Create a HDF5 group at a specified place in the file tree.
Definition: Hdf5File.cpp:178
Class wrapping the HDF5 routines.
Definition: Hdf5File.h:490
Hdf5File & operator=(const Hdf5File &)=delete
Operator = is not allowed.
void writeHyperSlab(const hid_t dataset, const DimensionSizes &position, const DimensionSizes &size, const T *data)
Write a hyperslab into the dataset.
Definition: Hdf5File.cpp:335
const std::string MatrixName
Datatype for matrix names.
Definition: MatrixNames.h:39
Structure with 4D dimension sizes (3 in space and 1 in time).
The matrix is stored in floating point 32b wide format.
The header file storing names of all variables.
MatrixDataType readMatrixDataType(const hid_t parentGroup, MatrixName &datasetName)
Read matrix data type from the dataset at a specified place in the file tree.
Definition: Hdf5File.cpp:821
The header file containing the structure with 3D dimension sizes.
The matrix is stored in fixed point point 64b wide format.
static const std::string kMatrixDataTypeNames[]
String representation of different data types.
Definition: Hdf5File.h:865
void closeGroup(const hid_t group)
Close a group.
Definition: Hdf5File.cpp:214
void close()
Close the HDF5 file.
Definition: Hdf5File.cpp:161
MatrixDomainType
HDF5 Matrix domain type (real or complex).
Definition: Hdf5File.h:512
void readCompleteDataset(const hid_t parentGroup, MatrixName &datasetName, const DimensionSizes &dimensionSizes, T *data)
Read data from the dataset at a specified place in the file tree.
Definition: Hdf5File.cpp:678
static const std::string kMatrixDataTypeName
String representation of the Data type in the HDF5 file.
Definition: Hdf5File.h:860
static const std::string kMatrixDomainTypeNames[]
String representation of different domain types.
Definition: Hdf5File.h:863
MatrixDataType
HDF5 matrix data type (float or uint64).
Definition: Hdf5File.h:499
Hdf5File()
Constructor of the class.
Definition: Hdf5File.cpp:73
hid_t createDataset(const hid_t parentGroup, MatrixName &datasetName, const DimensionSizes &dimensionSizes, const DimensionSizes &chunkSizes, const MatrixDataType matrixDataType, const size_t compressionLevel)
Create a float HDF5 dataset at a specified place in the file tree (3D/4D).
Definition: Hdf5File.cpp:241
bool isOpen() const
Is the file opened?
Definition: Hdf5File.h:560
MatrixDomainType readMatrixDomainType(const hid_t parentGroup, MatrixName &datasetName)
Read matrix dataset domain type at a specified place in the file tree.
Definition: Hdf5File.cpp:848