From d2a111dc40fe0dfeea7f9ce512421ed7b4badbb0 Mon Sep 17 00:00:00 2001 From: hinder Date: Wed, 21 Nov 2012 16:12:15 +0000 Subject: utils.cc: Add error checking to HDF5 calls git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinAnalysis/Multipole/trunk@87 4f5cb9a8-4dd8-4c2d-9bbd-173fa4467843 --- src/utils.cc | 41 ++++++++++++++++++++++++++++------------- 1 file changed, 28 insertions(+), 13 deletions(-) diff --git a/src/utils.cc b/src/utils.cc index 9b791b0..afd4278 100644 --- a/src/utils.cc +++ b/src/utils.cc @@ -25,6 +25,22 @@ #include "utils.hh" #include "integrate.hh" +// check return code of HDF5 call abort with an error message if there was an error. +// adapted from CarpetIOHDF5/src/CarpetIOHDF5.hh. +#define HDF5_ERROR(fn_call) \ + do { \ + int _error_code = fn_call; \ + \ + \ + if (_error_code < 0) \ + { \ + CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING, \ + "HDF5 call '%s' returned error code %d", \ + #fn_call, _error_code); \ + } \ + } while (0) + + using namespace std; //////////////////////////////////////////////////////////////////////// @@ -188,7 +204,6 @@ void Multipole_OutputComplexToH5File(CCTK_ARGUMENTS, const string &basename, con DECLARE_CCTK_PARAMETERS; string output_name = out_dir + string("/") + basename; - int status = 0; static map checked; // Has the given file been checked // for truncation? map<*,bool> // defaults to false @@ -221,7 +236,7 @@ void Multipole_OutputComplexToH5File(CCTK_ARGUMENTS, const string &basename, con hid_t cparms = -1; hsize_t chunk_dims[2] = {hdf5_chunk_size,3}; cparms = H5Pcreate (H5P_DATASET_CREATE); - status = H5Pset_chunk(cparms, 2, chunk_dims); + HDF5_ERROR(H5Pset_chunk(cparms, 2, chunk_dims)); dataset = H5Dcreate(file, datasetname.c_str(), H5T_NATIVE_DOUBLE, dataspace, cparms); H5Pclose(cparms); @@ -230,33 +245,33 @@ void Multipole_OutputComplexToH5File(CCTK_ARGUMENTS, const string &basename, con hid_t filespace = H5Dget_space (dataset); hsize_t filedims[2]; hsize_t maxdims[2]; - status = H5Sget_simple_extent_dims(filespace, filedims, maxdims); + HDF5_ERROR(H5Sget_simple_extent_dims(filespace, filedims, maxdims)); filedims[0] += 1; hsize_t size[2] = {filedims[0], filedims[1]}; - status = H5Dextend (dataset, size); - H5Sclose(filespace); + HDF5_ERROR(H5Dextend (dataset, size)); + HDF5_ERROR(H5Sclose(filespace)); /* Select a hyperslab */ hsize_t offset[2] = {filedims[0]-1, 0}; hsize_t dims2[2] = {1,3}; filespace = H5Dget_space (dataset); - status = H5Sselect_hyperslab (filespace, H5S_SELECT_SET, offset, NULL, - dims2, NULL); + HDF5_ERROR(H5Sselect_hyperslab (filespace, H5S_SELECT_SET, offset, NULL, + dims2, NULL)); CCTK_REAL data[] = {cctk_time, redata, imdata}; hid_t memdataspace = H5Screate_simple(2, dims2, NULL); /* Write the data to the hyperslab */ - status = H5Dwrite (dataset, H5T_NATIVE_DOUBLE, memdataspace, filespace, - H5P_DEFAULT, data); + HDF5_ERROR(H5Dwrite (dataset, H5T_NATIVE_DOUBLE, memdataspace, filespace, + H5P_DEFAULT, data)); - H5Dclose(dataset); - H5Sclose(filespace); - H5Sclose(memdataspace); + HDF5_ERROR(H5Dclose(dataset)); + HDF5_ERROR(H5Sclose(filespace)); + HDF5_ERROR(H5Sclose(memdataspace)); - status = H5Fclose(file); + HDF5_ERROR(H5Fclose(file)); } #else -- cgit v1.2.3