GADGET-4
hdf5_util.cc
Go to the documentation of this file.
1/*******************************************************************************
2 * \copyright This file is part of the GADGET4 N-body/SPH code developed
3 * \copyright by Volker Springel. Copyright (C) 2014-2020 by Volker Springel
4 * \copyright (vspringel@mpa-garching.mpg.de) and all contributing authors.
5 *******************************************************************************/
6
12#include "gadgetconfig.h"
13
14#include <hdf5.h>
15#include <math.h>
16#include <mpi.h>
17#include <stdio.h>
18#include <stdlib.h>
19#include <string.h>
20
21#include "../data/allvars.h"
22#include "../data/dtypes.h"
23#include "../io/hdf5_util.h"
24
25#define HALF_ROUND_STYLE 1
26#include "../half/half.hpp"
28
32
33static herr_t hdf5_conv_half(hid_t src, hid_t dst, H5T_cdata_t *cdata, size_t nelmts, size_t buf_str, size_t bkg_str, void *buf,
34 void *background, hid_t plist)
35{
36 size_t src_size = H5Tget_size(src);
37 size_t dst_size = H5Tget_size(dst);
38
39 char *src_buf = (char *)buf;
40 char *dst_buf = (char *)buf;
41
42 int direction;
43
44 switch(cdata->command)
45 {
46 case H5T_CONV_INIT:
47 /*
48 * We are being queried to see if we handle this
49 * conversion.
50 */
51 if(H5Tequal(src, Halfprec_memtype) || H5Tequal(dst, Halfprec_memtype))
52 {
53 cdata->need_bkg = H5T_BKG_NO;
54 return 0;
55 }
56 else
57 return -1;
58 break;
59
60 case H5T_CONV_FREE:
61 break;
62
63 case H5T_CONV_CONV:
64 /*
65 * Convert each element, watch out for overlap src
66 * with dst on the left-most element of the buffer.
67 * If the destination size is larger than the source size,
68 * then we must process the elements from right to left.
69 */
70
71 if(dst_size > src_size)
72 {
73 direction = -1;
74 src_buf += (nelmts - 1) * src_size;
75 dst_buf += (nelmts - 1) * dst_size;
76 }
77 else
78 {
79 direction = 1;
80 }
81
82 for(size_t i = 0; i < nelmts; i++)
83 {
84 if(src_size == 2)
85 {
86 if(dst_size == 4)
87 {
88 *((float *)dst_buf) = (float)(*((half *)src_buf));
89 }
90 else if(dst_size == 8)
91 {
92 *((double *)dst_buf) = (float)(*((half *)src_buf));
93 }
94 }
95 else if(dst_size == 2)
96 {
97 if(src_size == 4)
98 {
99 *((half *)dst_buf) = (half)(*((float *)src_buf));
100 }
101 else if(src_size == 8)
102 {
103 *((half *)dst_buf) = (half)(*((double *)src_buf));
104 }
105 }
106 src_buf += src_size * direction;
107 dst_buf += dst_size * direction;
108 }
109
110 break;
111
112 default:
113 /*
114 * Unknown command.
115 */
116 return -1;
117 }
118 return 0;
119}
120
122{
123 /* define the half-precision type */
124
125 Halfprec_memtype = H5Tcopy(H5T_NATIVE_FLOAT);
126
127 H5Tset_fields(Halfprec_memtype, 15, 10, 5, 0, 10);
128 H5Tset_ebias(Halfprec_memtype, 15);
129 H5Tset_precision(Halfprec_memtype, 16);
130 H5Tset_size(Halfprec_memtype, 2);
131
132 H5Tregister(H5T_PERS_SOFT, "half-converter", Halfprec_memtype, Halfprec_memtype, hdf5_conv_half);
133}
134
136{
137 Int48_memtype = H5Tcopy(H5T_NATIVE_UINT32);
138 H5Tset_precision(Int48_memtype, 48);
139
140 Int128_memtype = H5Tcopy(H5T_NATIVE_UINT64);
141 H5Tset_precision(Int128_memtype, 128);
142}
143
158hid_t my_H5Fcreate(const char *fname, unsigned int flags, hid_t fcpl_id, hid_t fapl_id)
159{
160 hid_t file_id = H5Fcreate(fname, flags, fcpl_id, fapl_id);
161
162 if(file_id < 0)
163 {
164 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
165 Terminate("Error detected in HDF5: unable to create file %s\n", fname);
166 }
167
168 return file_id;
169}
170
174hid_t my_H5Gcreate(hid_t loc_id, const char *groupname, size_t size_hint)
175{
176 hid_t group_id = H5Gcreate(loc_id, groupname, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
177
178 if(group_id < 0)
179 {
180 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
181 Terminate("Error detected in HDF5: unable to create group %s\n", groupname);
182 }
183
184 return group_id;
185}
186
190hid_t my_H5Dcreate(hid_t loc_id, const char *datasetname, hid_t type_id, hid_t space_id, hid_t dcpl_id)
191{
192 hid_t dataset_id = H5Dcreate(loc_id, datasetname, type_id, space_id, H5P_DEFAULT, dcpl_id, H5P_DEFAULT);
193
194 if(dataset_id < 0)
195 {
196 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
197 Terminate("Error detected in HDF5: unable to create dataset %s\n", datasetname);
198 }
199
200 return dataset_id;
201}
202
206herr_t my_H5Dwrite(hid_t dataset_id, hid_t mem_type_id, hid_t mem_space_id, hid_t file_space_id, hid_t xfer_plist_id, const void *buf,
207 const char *datasetname)
208{
209 herr_t status = H5Dwrite(dataset_id, mem_type_id, mem_space_id, file_space_id, xfer_plist_id, buf);
210
211 if(status < 0)
212 {
213 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
214 Terminate("Error detected in HDF5: unable to write dataset %s\n", datasetname);
215 }
216
217 return status;
218}
219
223hid_t my_H5Acreate(hid_t loc_id, const char *attr_name, hid_t type_id, hid_t space_id, hid_t acpl_id)
224{
225 if(H5Aexists(loc_id, attr_name))
226 H5Adelete(loc_id, attr_name);
227
228 hid_t attribute_id = H5Acreate(loc_id, attr_name, type_id, space_id, acpl_id, H5P_DEFAULT);
229
230 if(attribute_id < 0)
231 {
232 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
233 Terminate("Error detected in HDF5: unable to create attribute %s\n", attr_name);
234 }
235
236 return attribute_id;
237}
238
242herr_t my_H5Awrite(hid_t attr_id, hid_t mem_type_id, const void *buf, const char *attr_name)
243{
244 herr_t status = H5Awrite(attr_id, mem_type_id, buf);
245
246 if(status < 0)
247 {
248 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
249 Terminate("Error detected in HDF5: unable to write attribute %s\n", attr_name);
250 }
251
252 return status;
253}
254
258hid_t my_H5Screate(H5S_class_t type)
259{
260 hid_t dataspace_id = H5Screate(type);
261 if(dataspace_id < 0)
262 {
263 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
264 switch(type)
265 {
266 case H5S_SCALAR:
267 Terminate("Error detected in HDF5: unable to create a scalar dataspace\n");
268 break;
269 case H5S_SIMPLE:
270 Terminate("Error detected in HDF5: unable to create a simple dataspace\n");
271 break;
272 default:
273 Terminate("Error detected in HDF5: unknown dataspace type\n");
274 break;
275 }
276 }
277 return dataspace_id;
278}
279
283hid_t my_H5Screate_simple(int rank, const hsize_t *current_dims, const hsize_t *maximum_dims)
284{
285 hid_t dataspace_id = H5Screate_simple(rank, current_dims, maximum_dims);
286 if(dataspace_id < 0)
287 {
288 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
289 Terminate("Error detected in HDF5: unable to create a simple dataspace\n");
290 }
291 return dataspace_id;
292}
293
297hid_t my_H5Fopen(const char *fname, unsigned int flags, hid_t fapl_id)
298{
299 hid_t file_id = H5Fopen(fname, flags, fapl_id);
300 if(file_id < 0)
301 {
302 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
303 Terminate("Error detected in HDF5: unable to open file %s\n", fname);
304 }
305 return file_id;
306}
307
311hid_t my_H5Gopen(hid_t loc_id, const char *groupname)
312{
313 hid_t group = H5Gopen(loc_id, groupname, H5P_DEFAULT);
314 if(group < 0)
315 {
316 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
317 Terminate("Error detected in HDF5: unable to open group %s\n", groupname);
318 }
319 return group;
320}
321
325hid_t my_H5Dopen(hid_t file_id, const char *datasetname)
326{
327 hid_t dataset = H5Dopen(file_id, datasetname, H5P_DEFAULT);
328 if(dataset < 0)
329 {
330 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
331 Terminate("Error detected in HDF5: unable to open dataset %s\n", datasetname);
332 }
333 return dataset;
334}
335
336herr_t my_H5Dset_extent(hid_t dset_id, const hsize_t size[])
337{
338 herr_t status = H5Dset_extent(dset_id, size);
339 if(status < 0)
340 {
341 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
342 Terminate("Error detected in HDF5: unable to set extent of dataset\n");
343 }
344 return status;
345}
346
352hid_t my_H5Dopen_if_existing(hid_t file_id, const char *datasetname)
353{
354 /* save error handler and disable it */
355 H5E_auto_t errfunc;
356 void *client_data;
357 H5Eget_auto(H5E_DEFAULT, &errfunc, &client_data);
358 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
359
360 hid_t dataset = H5Dopen(file_id, datasetname, H5P_DEFAULT);
361
362 /* reset error handler */
363 H5Eset_auto(H5E_DEFAULT, errfunc, client_data);
364
365 return dataset;
366}
367
371hid_t my_H5Aopen_name(hid_t loc_id, const char *attr_name)
372{
373 hid_t attribute_id = H5Aopen_name(loc_id, attr_name);
374 if(attribute_id < 0)
375 {
376 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
377 Terminate("Error detected in HDF5: unable to open attribute %s\n", attr_name);
378 }
379 return attribute_id;
380}
381
385herr_t my_H5Dread(hid_t dataset_id, hid_t mem_type_id, hid_t mem_space_id, hid_t file_space_id, hid_t xfer_plist_id, void *buf,
386 const char *datasetname)
387{
388 herr_t status = H5Dread(dataset_id, mem_type_id, mem_space_id, file_space_id, xfer_plist_id, buf);
389 if(status < 0)
390 {
391 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
392 Terminate("Error detected in HDF5: unable to read dataset %s\n", datasetname);
393 }
394 return status;
395}
396
400herr_t my_H5Aread(hid_t attr_id, hid_t mem_type_id, void *buf, const char *attr_name, hssize_t size)
401{
402 hid_t hdf5_space = H5Aget_space(attr_id);
403 hssize_t attr_size = H5Sget_simple_extent_npoints(hdf5_space);
404 H5Sclose(hdf5_space);
405
406 if(attr_size != size)
407 {
408 H5E_auto_t errfunc;
409 void *client_data;
410 H5Eget_auto(H5E_DEFAULT, &errfunc, &client_data);
411 errfunc(H5P_DEFAULT, client_data);
412 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
413 Terminate("Error detected in HDF5: mismatch in size for attribute %s, expected size = %lld, actual attribute size = %lld\n",
414 attr_name, size, attr_size);
415 }
416
417 herr_t status = H5Aread(attr_id, mem_type_id, buf);
418 if(status < 0)
419 {
420 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
421 Terminate("Error detected in HDF5: unable to read attribute %s\n", attr_name);
422 }
423 return status;
424}
425
429herr_t my_H5Aclose(hid_t attr_id, const char *attr_name)
430{
431 herr_t status = H5Aclose(attr_id);
432 if(status < 0)
433 {
434 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
435 Terminate("Error detected in HDF5: unable to close attribute %s\n", attr_name);
436 }
437 return status;
438}
439
443herr_t my_H5Dclose(hid_t dataset_id, const char *datasetname)
444{
445 herr_t status = H5Dclose(dataset_id);
446 if(status < 0)
447 {
448 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
449 Terminate("Error detected in HDF5: unable to close dataset %s\n", datasetname);
450 }
451 return status;
452}
453
457herr_t my_H5Gclose(hid_t group_id, const char *groupname)
458{
459 herr_t status = H5Gclose(group_id);
460 if(status < 0)
461 {
462 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
463 Terminate("Error detected in HDF5: unable to close group %s\n", groupname);
464 }
465 return status;
466}
467
468herr_t my_H5Pclose(hid_t plist)
469{
470 herr_t status = H5Pclose(plist);
471 if(status < 0)
472 {
473 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
474 Terminate("Error detected in HDF5: unable to close property\n");
475 }
476 return status;
477}
478
482herr_t my_H5Fclose(hid_t file_id, const char *fname)
483{
484 herr_t status = H5Fclose(file_id);
485 if(status < 0)
486 {
487 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
488 Terminate("Error detected in HDF5: unable to close file %s\n", fname);
489 }
490 return status;
491}
492
496herr_t my_H5Sclose(hid_t dataspace_id, H5S_class_t type)
497{
498 herr_t status = H5Sclose(dataspace_id);
499 if(status < 0)
500 {
501 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
502 switch(type)
503 {
504 case H5S_SCALAR:
505 Terminate("Error detected in HDF5: unable to close a scalar dataspace\n");
506 break;
507 case H5S_SIMPLE:
508 Terminate("Error detected in HDF5: unable to close a simple dataspace\n");
509 break;
510 default:
511 Terminate("Error detected in HDF5: unknown dataspace type\n");
512 break;
513 }
514 }
515 return status;
516}
517
521hid_t my_H5Tcopy(hid_t type_id)
522{
523 hid_t datatype_id = H5Tcopy(type_id);
524 if(datatype_id < 0)
525 {
526 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
527 Terminate("Error detected in HDF5: could not properly copy datatype\n");
528 }
529 return datatype_id;
530}
531
535herr_t my_H5Tclose(hid_t type_id)
536{
537 herr_t status = H5Tclose(type_id);
538 if(status < 0)
539 {
540 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
541 Terminate("Error detected in HDF5: could not properly close datatype\n");
542 }
543 return status;
544}
545
549herr_t my_H5Sselect_hyperslab(hid_t space_id, H5S_seloper_t op, const hsize_t *start, const hsize_t *stride, const hsize_t *count,
550 const hsize_t *block)
551{
552 herr_t status = H5Sselect_hyperslab(space_id, op, start, stride, count, block);
553 if(status < 0)
554 {
555 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
556 Terminate("Error detected in HDF5: could not properly select the chosen hyperslab\n");
557 }
558 return status;
559}
560
564size_t my_H5Tget_size(hid_t datatype_id)
565{
566 size_t size = H5Tget_size(datatype_id);
567 if(size == 0)
568 {
569 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
570 Terminate("Error detected in HDF5: unable to determine the size of the given datatype\n");
571 }
572 return size;
573}
574
578herr_t my_H5Tset_size(hid_t datatype_id, size_t size)
579{
580 herr_t status = H5Tset_size(datatype_id, size);
581 if(status < 0)
582 {
583 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
584 Terminate("Error detected in HDF5: could not properly set the size of the given datatype\n");
585 }
586 return status;
587}
588
589void read_scalar_attribute(hid_t handle, const char *attr_name, void *buf, hid_t mem_type_id)
590{
591 hid_t hdf5_attribute = my_H5Aopen_name(handle, attr_name);
592 my_H5Aread(hdf5_attribute, mem_type_id, buf, attr_name, 1);
593 my_H5Aclose(hdf5_attribute, attr_name);
594}
595
596int read_scalar_attribute(hid_t handle, const char *attr_name, const char *alternative_name, void *buf, hid_t mem_type_id)
597{
598 if(H5Aexists(handle, attr_name))
599 {
600 hid_t hdf5_attribute = my_H5Aopen_name(handle, attr_name);
601 my_H5Aread(hdf5_attribute, mem_type_id, buf, attr_name, 1);
602 my_H5Aclose(hdf5_attribute, attr_name);
603 return 0;
604 }
605 else
606 {
607 hid_t hdf5_attribute = my_H5Aopen_name(handle, alternative_name);
608 my_H5Aread(hdf5_attribute, mem_type_id, buf, alternative_name, 1);
609 my_H5Aclose(hdf5_attribute, alternative_name);
610 return 1;
611 }
612}
613
614void read_vector_attribute(hid_t handle, const char *attr_name, void *buf, hid_t mem_type_id, int length)
615{
616 hid_t hdf5_attribute = my_H5Aopen_name(handle, attr_name);
617 my_H5Aread(hdf5_attribute, mem_type_id, buf, attr_name, length);
618 my_H5Aclose(hdf5_attribute, attr_name);
619}
620
621void write_scalar_attribute(hid_t handle, const char *attr_name, const void *buf, hid_t mem_type_id)
622{
623 if(H5Aexists(handle, attr_name))
624 H5Adelete(handle, attr_name);
625
626 hid_t hdf5_dataspace = my_H5Screate(H5S_SCALAR);
627
628 hid_t hdf5_attribute = my_H5Acreate(handle, attr_name, mem_type_id, hdf5_dataspace, H5P_DEFAULT);
629
630 my_H5Awrite(hdf5_attribute, mem_type_id, buf, attr_name);
631
632 my_H5Aclose(hdf5_attribute, attr_name);
633 my_H5Sclose(hdf5_dataspace, H5S_SCALAR);
634}
635
636void write_vector_attribute(hid_t handle, const char *attr_name, const void *buf, hid_t mem_type_id, int length)
637{
638 if(H5Aexists(handle, attr_name))
639 H5Adelete(handle, attr_name);
640
641 hsize_t adim[1] = {(hsize_t)length};
642
643 hid_t hdf5_dataspace = my_H5Screate(H5S_SIMPLE);
644 H5Sset_extent_simple(hdf5_dataspace, 1, adim, NULL);
645
646 hid_t hdf5_attribute = my_H5Acreate(handle, attr_name, mem_type_id, hdf5_dataspace, H5P_DEFAULT);
647
648 my_H5Awrite(hdf5_attribute, mem_type_id, buf, attr_name);
649
650 my_H5Aclose(hdf5_attribute, attr_name);
651 my_H5Sclose(hdf5_dataspace, H5S_SIMPLE);
652}
653
654void write_string_attribute(hid_t handle, const char *attr_name, const char *buf)
655{
656 if(H5Aexists(handle, attr_name))
657 H5Adelete(handle, attr_name);
658
659 hid_t atype = my_H5Tcopy(H5T_C_S1);
660 my_H5Tset_size(atype, strlen(buf));
661
662 hid_t hdf5_dataspace = my_H5Screate(H5S_SCALAR);
663 hid_t hdf5_attribute = my_H5Acreate(handle, attr_name, atype, hdf5_dataspace, H5P_DEFAULT);
664
665 my_H5Awrite(hdf5_attribute, atype, buf, attr_name);
666
667 my_H5Aclose(hdf5_attribute, attr_name);
668 my_H5Sclose(hdf5_dataspace, H5S_SCALAR);
669}
void write_vector_attribute(hid_t handle, const char *attr_name, const void *buf, hid_t mem_type_id, int length)
Definition: hdf5_util.cc:636
hid_t my_H5Gopen(hid_t loc_id, const char *groupname)
Definition: hdf5_util.cc:311
herr_t my_H5Sselect_hyperslab(hid_t space_id, H5S_seloper_t op, const hsize_t *start, const hsize_t *stride, const hsize_t *count, const hsize_t *block)
Definition: hdf5_util.cc:549
hid_t my_H5Dopen(hid_t file_id, const char *datasetname)
Definition: hdf5_util.cc:325
herr_t my_H5Dread(hid_t dataset_id, hid_t mem_type_id, hid_t mem_space_id, hid_t file_space_id, hid_t xfer_plist_id, void *buf, const char *datasetname)
Definition: hdf5_util.cc:385
hid_t my_H5Fopen(const char *fname, unsigned int flags, hid_t fapl_id)
Definition: hdf5_util.cc:297
void my_create_HDF5_special_integer_types(void)
Definition: hdf5_util.cc:135
herr_t my_H5Gclose(hid_t group_id, const char *groupname)
Definition: hdf5_util.cc:457
hid_t Halfprec_memtype
Definition: hdf5_util.cc:29
size_t my_H5Tget_size(hid_t datatype_id)
Definition: hdf5_util.cc:564
herr_t my_H5Tclose(hid_t type_id)
Definition: hdf5_util.cc:535
herr_t my_H5Tset_size(hid_t datatype_id, size_t size)
Definition: hdf5_util.cc:578
herr_t my_H5Dset_extent(hid_t dset_id, const hsize_t size[])
Definition: hdf5_util.cc:336
void read_vector_attribute(hid_t handle, const char *attr_name, void *buf, hid_t mem_type_id, int length)
Definition: hdf5_util.cc:614
hid_t Int48_memtype
Definition: hdf5_util.cc:30
hid_t my_H5Gcreate(hid_t loc_id, const char *groupname, size_t size_hint)
Definition: hdf5_util.cc:174
herr_t my_H5Pclose(hid_t plist)
Definition: hdf5_util.cc:468
hid_t my_H5Screate_simple(int rank, const hsize_t *current_dims, const hsize_t *maximum_dims)
Definition: hdf5_util.cc:283
herr_t my_H5Fclose(hid_t file_id, const char *fname)
Definition: hdf5_util.cc:482
herr_t my_H5Dclose(hid_t dataset_id, const char *datasetname)
Definition: hdf5_util.cc:443
hid_t my_H5Fcreate(const char *fname, unsigned int flags, hid_t fcpl_id, hid_t fapl_id)
Definition: hdf5_util.cc:158
hid_t my_H5Aopen_name(hid_t loc_id, const char *attr_name)
Definition: hdf5_util.cc:371
hid_t my_H5Dopen_if_existing(hid_t file_id, const char *datasetname)
Definition: hdf5_util.cc:352
herr_t my_H5Aread(hid_t attr_id, hid_t mem_type_id, void *buf, const char *attr_name, hssize_t size)
Definition: hdf5_util.cc:400
herr_t my_H5Aclose(hid_t attr_id, const char *attr_name)
Definition: hdf5_util.cc:429
hid_t my_H5Dcreate(hid_t loc_id, const char *datasetname, hid_t type_id, hid_t space_id, hid_t dcpl_id)
Definition: hdf5_util.cc:190
herr_t my_H5Sclose(hid_t dataspace_id, H5S_class_t type)
Definition: hdf5_util.cc:496
hid_t Int128_memtype
Definition: hdf5_util.cc:31
herr_t my_H5Awrite(hid_t attr_id, hid_t mem_type_id, const void *buf, const char *attr_name)
Definition: hdf5_util.cc:242
void write_scalar_attribute(hid_t handle, const char *attr_name, const void *buf, hid_t mem_type_id)
Definition: hdf5_util.cc:621
hid_t my_H5Screate(H5S_class_t type)
Definition: hdf5_util.cc:258
herr_t my_H5Dwrite(hid_t dataset_id, hid_t mem_type_id, hid_t mem_space_id, hid_t file_space_id, hid_t xfer_plist_id, const void *buf, const char *datasetname)
Definition: hdf5_util.cc:206
hid_t my_H5Acreate(hid_t loc_id, const char *attr_name, hid_t type_id, hid_t space_id, hid_t acpl_id)
Definition: hdf5_util.cc:223
void read_scalar_attribute(hid_t handle, const char *attr_name, void *buf, hid_t mem_type_id)
Definition: hdf5_util.cc:589
hid_t my_H5Tcopy(hid_t type_id)
Definition: hdf5_util.cc:521
void my_create_HDF5_halfprec_handler(void)
Definition: hdf5_util.cc:121
void write_string_attribute(hid_t handle, const char *attr_name, const char *buf)
Definition: hdf5_util.cc:654
#define Terminate(...)
Definition: macros.h:15