12#include "gadgetconfig.h"
21#include "../data/allvars.h"
22#include "../data/dtypes.h"
23#include "../io/hdf5_util.h"
25#define HALF_ROUND_STYLE 1
26#include "../half/half.hpp"
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)
36 size_t src_size = H5Tget_size(src);
37 size_t dst_size = H5Tget_size(dst);
39 char *src_buf = (
char *)buf;
40 char *dst_buf = (
char *)buf;
44 switch(cdata->command)
53 cdata->need_bkg = H5T_BKG_NO;
71 if(dst_size > src_size)
74 src_buf += (nelmts - 1) * src_size;
75 dst_buf += (nelmts - 1) * dst_size;
82 for(
size_t i = 0; i < nelmts; i++)
88 *((
float *)dst_buf) = (float)(*((
half *)src_buf));
90 else if(dst_size == 8)
92 *((
double *)dst_buf) = (float)(*((
half *)src_buf));
95 else if(dst_size == 2)
99 *((
half *)dst_buf) = (
half)(*((
float *)src_buf));
101 else if(src_size == 8)
103 *((
half *)dst_buf) = (
half)(*((
double *)src_buf));
106 src_buf += src_size * direction;
107 dst_buf += dst_size * direction;
158hid_t
my_H5Fcreate(
const char *fname,
unsigned int flags, hid_t fcpl_id, hid_t fapl_id)
160 hid_t file_id = H5Fcreate(fname, flags, fcpl_id, fapl_id);
164 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
165 Terminate(
"Error detected in HDF5: unable to create file %s\n", fname);
174hid_t
my_H5Gcreate(hid_t loc_id,
const char *groupname,
size_t size_hint)
176 hid_t group_id = H5Gcreate(loc_id, groupname, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
180 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
181 Terminate(
"Error detected in HDF5: unable to create group %s\n", groupname);
190hid_t
my_H5Dcreate(hid_t loc_id,
const char *datasetname, hid_t type_id, hid_t space_id, hid_t dcpl_id)
192 hid_t dataset_id = H5Dcreate(loc_id, datasetname, type_id, space_id, H5P_DEFAULT, dcpl_id, H5P_DEFAULT);
196 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
197 Terminate(
"Error detected in HDF5: unable to create dataset %s\n", datasetname);
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)
209 herr_t status = H5Dwrite(dataset_id, mem_type_id, mem_space_id, file_space_id, xfer_plist_id, buf);
213 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
214 Terminate(
"Error detected in HDF5: unable to write dataset %s\n", datasetname);
223hid_t
my_H5Acreate(hid_t loc_id,
const char *attr_name, hid_t type_id, hid_t space_id, hid_t acpl_id)
225 if(H5Aexists(loc_id, attr_name))
226 H5Adelete(loc_id, attr_name);
228 hid_t attribute_id = H5Acreate(loc_id, attr_name, type_id, space_id, acpl_id, H5P_DEFAULT);
232 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
233 Terminate(
"Error detected in HDF5: unable to create attribute %s\n", attr_name);
242herr_t
my_H5Awrite(hid_t attr_id, hid_t mem_type_id,
const void *buf,
const char *attr_name)
244 herr_t status = H5Awrite(attr_id, mem_type_id, buf);
248 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
249 Terminate(
"Error detected in HDF5: unable to write attribute %s\n", attr_name);
260 hid_t dataspace_id = H5Screate(type);
263 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
267 Terminate(
"Error detected in HDF5: unable to create a scalar dataspace\n");
270 Terminate(
"Error detected in HDF5: unable to create a simple dataspace\n");
273 Terminate(
"Error detected in HDF5: unknown dataspace type\n");
285 hid_t dataspace_id = H5Screate_simple(rank, current_dims, maximum_dims);
288 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
289 Terminate(
"Error detected in HDF5: unable to create a simple dataspace\n");
297hid_t
my_H5Fopen(
const char *fname,
unsigned int flags, hid_t fapl_id)
299 hid_t file_id = H5Fopen(fname, flags, fapl_id);
302 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
303 Terminate(
"Error detected in HDF5: unable to open file %s\n", fname);
313 hid_t group = H5Gopen(loc_id, groupname, H5P_DEFAULT);
316 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
317 Terminate(
"Error detected in HDF5: unable to open group %s\n", groupname);
327 hid_t dataset = H5Dopen(file_id, datasetname, H5P_DEFAULT);
330 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
331 Terminate(
"Error detected in HDF5: unable to open dataset %s\n", datasetname);
338 herr_t status = H5Dset_extent(dset_id, size);
341 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
342 Terminate(
"Error detected in HDF5: unable to set extent of dataset\n");
357 H5Eget_auto(H5E_DEFAULT, &errfunc, &client_data);
358 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
360 hid_t dataset = H5Dopen(file_id, datasetname, H5P_DEFAULT);
363 H5Eset_auto(H5E_DEFAULT, errfunc, client_data);
373 hid_t attribute_id = H5Aopen_name(loc_id, attr_name);
376 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
377 Terminate(
"Error detected in HDF5: unable to open attribute %s\n", attr_name);
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)
388 herr_t status = H5Dread(dataset_id, mem_type_id, mem_space_id, file_space_id, xfer_plist_id, buf);
391 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
392 Terminate(
"Error detected in HDF5: unable to read dataset %s\n", datasetname);
400herr_t
my_H5Aread(hid_t attr_id, hid_t mem_type_id,
void *buf,
const char *attr_name, hssize_t size)
402 hid_t hdf5_space = H5Aget_space(attr_id);
403 hssize_t attr_size = H5Sget_simple_extent_npoints(hdf5_space);
404 H5Sclose(hdf5_space);
406 if(attr_size != size)
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);
417 herr_t status = H5Aread(attr_id, mem_type_id, buf);
420 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
421 Terminate(
"Error detected in HDF5: unable to read attribute %s\n", attr_name);
431 herr_t status = H5Aclose(attr_id);
434 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
435 Terminate(
"Error detected in HDF5: unable to close attribute %s\n", attr_name);
445 herr_t status = H5Dclose(dataset_id);
448 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
449 Terminate(
"Error detected in HDF5: unable to close dataset %s\n", datasetname);
459 herr_t status = H5Gclose(group_id);
462 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
463 Terminate(
"Error detected in HDF5: unable to close group %s\n", groupname);
470 herr_t status = H5Pclose(plist);
473 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
474 Terminate(
"Error detected in HDF5: unable to close property\n");
484 herr_t status = H5Fclose(file_id);
487 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
488 Terminate(
"Error detected in HDF5: unable to close file %s\n", fname);
498 herr_t status = H5Sclose(dataspace_id);
501 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
505 Terminate(
"Error detected in HDF5: unable to close a scalar dataspace\n");
508 Terminate(
"Error detected in HDF5: unable to close a simple dataspace\n");
511 Terminate(
"Error detected in HDF5: unknown dataspace type\n");
523 hid_t datatype_id = H5Tcopy(type_id);
526 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
527 Terminate(
"Error detected in HDF5: could not properly copy datatype\n");
537 herr_t status = H5Tclose(type_id);
540 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
541 Terminate(
"Error detected in HDF5: could not properly close datatype\n");
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)
552 herr_t status = H5Sselect_hyperslab(space_id, op, start, stride, count, block);
555 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
556 Terminate(
"Error detected in HDF5: could not properly select the chosen hyperslab\n");
566 size_t size = H5Tget_size(datatype_id);
569 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
570 Terminate(
"Error detected in HDF5: unable to determine the size of the given datatype\n");
580 herr_t status = H5Tset_size(datatype_id, size);
583 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
584 Terminate(
"Error detected in HDF5: could not properly set the size of the given datatype\n");
592 my_H5Aread(hdf5_attribute, mem_type_id, buf, attr_name, 1);
596int read_scalar_attribute(hid_t handle,
const char *attr_name,
const char *alternative_name,
void *buf, hid_t mem_type_id)
598 if(H5Aexists(handle, attr_name))
601 my_H5Aread(hdf5_attribute, mem_type_id, buf, attr_name, 1);
608 my_H5Aread(hdf5_attribute, mem_type_id, buf, alternative_name, 1);
617 my_H5Aread(hdf5_attribute, mem_type_id, buf, attr_name, length);
623 if(H5Aexists(handle, attr_name))
624 H5Adelete(handle, attr_name);
628 hid_t hdf5_attribute =
my_H5Acreate(handle, attr_name, mem_type_id, hdf5_dataspace, H5P_DEFAULT);
630 my_H5Awrite(hdf5_attribute, mem_type_id, buf, attr_name);
638 if(H5Aexists(handle, attr_name))
639 H5Adelete(handle, attr_name);
641 hsize_t adim[1] = {(hsize_t)length};
644 H5Sset_extent_simple(hdf5_dataspace, 1, adim, NULL);
646 hid_t hdf5_attribute =
my_H5Acreate(handle, attr_name, mem_type_id, hdf5_dataspace, H5P_DEFAULT);
648 my_H5Awrite(hdf5_attribute, mem_type_id, buf, attr_name);
656 if(H5Aexists(handle, attr_name))
657 H5Adelete(handle, attr_name);
663 hid_t hdf5_attribute =
my_H5Acreate(handle, attr_name, atype, hdf5_dataspace, H5P_DEFAULT);
665 my_H5Awrite(hdf5_attribute, atype, buf, attr_name);
void write_vector_attribute(hid_t handle, const char *attr_name, const void *buf, hid_t mem_type_id, int length)
hid_t my_H5Gopen(hid_t loc_id, const char *groupname)
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)
hid_t my_H5Dopen(hid_t file_id, const char *datasetname)
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)
hid_t my_H5Fopen(const char *fname, unsigned int flags, hid_t fapl_id)
void my_create_HDF5_special_integer_types(void)
herr_t my_H5Gclose(hid_t group_id, const char *groupname)
size_t my_H5Tget_size(hid_t datatype_id)
herr_t my_H5Tclose(hid_t type_id)
herr_t my_H5Tset_size(hid_t datatype_id, size_t size)
herr_t my_H5Dset_extent(hid_t dset_id, const hsize_t size[])
void read_vector_attribute(hid_t handle, const char *attr_name, void *buf, hid_t mem_type_id, int length)
hid_t my_H5Gcreate(hid_t loc_id, const char *groupname, size_t size_hint)
herr_t my_H5Pclose(hid_t plist)
hid_t my_H5Screate_simple(int rank, const hsize_t *current_dims, const hsize_t *maximum_dims)
herr_t my_H5Fclose(hid_t file_id, const char *fname)
herr_t my_H5Dclose(hid_t dataset_id, const char *datasetname)
hid_t my_H5Fcreate(const char *fname, unsigned int flags, hid_t fcpl_id, hid_t fapl_id)
hid_t my_H5Aopen_name(hid_t loc_id, const char *attr_name)
hid_t my_H5Dopen_if_existing(hid_t file_id, const char *datasetname)
herr_t my_H5Aread(hid_t attr_id, hid_t mem_type_id, void *buf, const char *attr_name, hssize_t size)
herr_t my_H5Aclose(hid_t attr_id, const char *attr_name)
hid_t my_H5Dcreate(hid_t loc_id, const char *datasetname, hid_t type_id, hid_t space_id, hid_t dcpl_id)
herr_t my_H5Sclose(hid_t dataspace_id, H5S_class_t type)
herr_t my_H5Awrite(hid_t attr_id, hid_t mem_type_id, const void *buf, const char *attr_name)
void write_scalar_attribute(hid_t handle, const char *attr_name, const void *buf, hid_t mem_type_id)
hid_t my_H5Screate(H5S_class_t type)
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)
hid_t my_H5Acreate(hid_t loc_id, const char *attr_name, hid_t type_id, hid_t space_id, hid_t acpl_id)
void read_scalar_attribute(hid_t handle, const char *attr_name, void *buf, hid_t mem_type_id)
hid_t my_H5Tcopy(hid_t type_id)
void my_create_HDF5_halfprec_handler(void)
void write_string_attribute(hid_t handle, const char *attr_name, const char *buf)