23 class HDF5File :
public std::enable_shared_from_this<HDF5File> {
33 HDF5File(std::string
const& filename_);
44 void Open(std::string
const& filename_);
56 void Copy(std::string
const& destName, std::shared_ptr<HDF5File> srcFile, std::string
const& srcName);
69 template<
typename scalarType,
int fixedRows,
int fixedCols>
70 void WriteMatrix(std::string name, Eigen::Matrix<scalarType, fixedRows, fixedCols>
const& dataset) {
73 if( name.at(0)!=
'/' ) {
74 std::cerr << std::endl <<
"ERROR: Paths in the HDF5 file must start with a forward slash (/)" << std::endl <<
"\tHDF5File::WriteMatrix(std::string const&, Eigen::Matrix<scalarType, fixedRows, fixedCols> const&)" << std::endl << std::endl;
75 assert(name.at(0) ==
'/');
78 if((name.at(0)==
'/')&&(name.size()>0)){
80 name.erase(name.begin());
88 const hsize_t maxdim[2] = {H5S_UNLIMITED, H5S_UNLIMITED};
90 const hsize_t dimsf[2] = {(hsize_t)dataset.rows(), (hsize_t)dataset.cols()};
93 const hid_t filespace = H5Screate_simple(2, dimsf, maxdim);
101 dataID = H5Dopen(
fileID, name.c_str(), H5P_DEFAULT);
105 H5Sget_simple_extent_dims(H5Dget_space(dataID), dims,
nullptr);
107 if( (dimsf[0]!=dims[0]) || (dimsf[1]!=dims[1]) ) {
109 H5Dset_extent(dataID, dimsf);
121 const hid_t prop = H5Pcreate(H5P_DATASET_CREATE);
124 hsize_t chunk[2] = {100,100};
126 chunk[0] = std::min<hsize_t>(chunk[0],fixedRows);
130 chunk[1] = std::min<hsize_t>(chunk[1],fixedRows);
132 H5Pset_chunk(prop, 2, chunk);
135 dataID = H5Dcreate(
fileID,
148 Eigen::Matrix<scalarType, fixedRows, fixedCols> dataset_ = dataset.transpose();
170 template<
typename scalarType,
int fixedRows,
int fixedCols>
172 Eigen::Matrix<scalarType, fixedRows, fixedCols>
const& data,
173 unsigned int const row,
unsigned int const col)
176 if( name.at(0)!=
'/' ) {
177 std::cerr << std::endl <<
"ERROR: Paths in the HDF5 file must start with a forward slash (/)" << std::endl <<
"\tHDF5File::WriteMatrix(std::string const&, Eigen::Matrix<scalarType, fixedRows, fixedCols> const&)" << std::endl << std::endl;
178 assert(name.at(0) ==
'/');
185 const hsize_t stride[2] = {1, 1};
188 const hsize_t count[2] = {1, 1};
192 std::cerr << std::endl <<
"ERROR: Dataset " << name <<
" does not exsts." << std::endl << std::endl;
197 const hid_t dataID = H5Dopen(
fileID, name.c_str(), H5P_DEFAULT);
200 const hid_t fullspace = H5Dget_space(dataID);
203 const hsize_t start[2] = {row, col};
206 const hsize_t block[2] = {(hsize_t)data.rows(), (hsize_t)data.cols()};
209 H5Sselect_hyperslab(fullspace, H5S_SELECT_SET, start, stride, count, block);
212 const hsize_t limitedBlock = H5Screate_simple(2, block, block);
215 const Eigen::Matrix<scalarType, fixedRows, fixedCols>& data_ = data.transpose();
218 auto write = H5Dwrite(dataID,
226 H5Sclose(limitedBlock);
238 template<
typename scalarType =
double,
int fixedRows = Eigen::Dynamic,
int fixedCols = Eigen::Dynamic>
239 Eigen::Matrix<scalarType, fixedRows, fixedCols>
ReadMatrix(std::string
const& name)
const {
246 return Eigen::Matrix<scalarType, fixedRows, fixedCols>();
250 const hid_t dataset = H5Dopen2(hid_t(
fileID), name.c_str(), H5P_DEFAULT);
253 const hid_t filespace = H5Dget_space(dataset);
255 const int ndims = H5Sget_simple_extent_ndims(filespace);
258 std::vector<hsize_t> dimsr(ndims);
259 H5Sget_simple_extent_dims(filespace, dimsr.data(),
nullptr);
262 Eigen::Matrix<scalarType, Eigen::Dynamic, Eigen::Dynamic> pt;
265 pt.resize(dimsr[1], dimsr[0]);
268 pt = Eigen::Matrix<scalarType, 1, fixedCols>(dimsr[0]);
270 pt = Eigen::Matrix<scalarType, fixedRows, 1>(dimsr[0]);
276 const hid_t prop = H5Dget_create_plist(dataset);
279 const hid_t memspace = H5Screate_simple(ndims, dimsr.data(),
nullptr);
301 return pt.transpose();
310 template<
typename scalarType=
double,
int fixedRows = Eigen::Dynamic,
int fixedCols = Eigen::Dynamic>
324 return Eigen::Matrix<scalarType, fixedRows, fixedCols>();
328 const hid_t dataset = H5Dopen2(
fileID, name.c_str(), H5P_DEFAULT);
331 const hid_t dataspace = H5Dget_space(dataset);
335 H5Sget_simple_extent_dims(dataspace, dimsr,
nullptr);
338 Eigen::Matrix<scalarType, fixedCols, fixedRows> pt;
340 pt.resize(numCols,numRows);
342 assert(numRows<=dimsr[0]);
343 assert(numCols<=dimsr[1]);
345 assert(numRows<dimsr[0]);
348 const hsize_t start[2] = {(hsize_t) rowStart, (hsize_t) colStart};
349 const hsize_t count[2] = {(hsize_t) numRows, (hsize_t) numCols};
351 herr_t status = H5Sselect_hyperslab(dataspace, H5S_SELECT_SET, start, NULL, count, NULL);
355 hid_t memspace = H5Screate_simple(2, dimsr,
nullptr);
371 return pt.transpose();
393 bool IsDataSet(std::string
const& name)
const;
394 bool IsGroup(std::string
const& name)
const;
418 template<
typename scalarType>
426 const hsize_t maxdim[2] = {H5S_UNLIMITED, H5S_UNLIMITED};
429 const hsize_t dimsf[2] = {(hsize_t)rows, (hsize_t)cols};
432 const hid_t filespace = H5Screate_simple(2, dimsf, maxdim);
441 dataID = H5Dopen(
fileID, name.c_str(), H5P_DEFAULT);
445 H5Sget_simple_extent_dims(H5Dget_space(dataID), dims,
nullptr);
447 if( (dimsf[0]!=dims[0]) || (dimsf[1]!=dims[1]) ) {
449 H5Dset_extent(dataID, dimsf);
461 const hid_t prop = H5Pcreate(H5P_DATASET_CREATE);
464 hsize_t chunk[2] = {std::min<hsize_t>(100,rows), std::min<hsize_t>(100,cols)};
465 H5Pset_chunk(prop, 2, chunk);
490 template<
typename scalarType,
int fixedRows>
492 std::string
const& attributeName,
493 Eigen::Matrix<scalarType, fixedRows, 1>
const& attribute)
506 attributeName.c_str(),
508 (
void*) attribute.data(),
510 assert(set_result>=0);
521 template<
typename attType =
double>
523 std::string
const& attributeName,
524 attType
const& attribute)
527 Eigen::Matrix<attType, 1, 1> vec(1);
531 WriteVectorAttribute<attType, 1>(datasetName, attributeName, vec);
543 std::string
const& attributeName,
544 std::string
const& attribute);
551 template<
typename scalarType>
553 std::string
const& attributeName)
const {
556 Eigen::Matrix<scalarType, Eigen::Dynamic, 1> att = GetVectorAttribute<scalarType>(datasetName, attributeName);
557 assert(att.size() == 1);
568 template<
typename scalarType>
570 std::string
const& attributeName)
const
581 auto get_ndim = H5LTget_attribute_ndims(
fileID, datasetName.c_str(), attributeName.c_str(), rank);
582 assert(get_ndim >= 0);
585 H5T_class_t typeClass;
587 auto status = H5LTget_attribute_info(
fileID, datasetName.c_str(), attributeName.c_str(), dims, &typeClass, &typeSize);
591 Eigen::Matrix<scalarType, Eigen::Dynamic, 1> data(dims[0]);
596 assert(get_att >= 0);
608 std::string
GetStringAttribute(std::string
const& datasetName, std::string
const& attributeName)
const;
612 std::vector<std::string>
GetChildren(std::string base =
"/")
const;
622 void MergeFile(std::shared_ptr<HDF5File>
const& otherFile);
A wrapper around HDF5 (and PHDF5) to handle file I/O.
Eigen::Matrix< scalarType, Eigen::Dynamic, 1 > GetVectorAttribute(std::string const &datasetName, std::string const &attributeName) const
Read a vector attribute from the HDF5 file.
void Open(std::string const &filename_)
Opens or creates the file.
HDF5File(std::string const &filename_)
Open or create the file.
bool DoesGroupExist(std::string const &name) const
Check to see if a group exists.
virtual ~HDF5File()
If HDF5File is destroyed, the file should be closed.
void Close()
Close the file.
hid_t fileID
The HDF5 file ID.
void CreateDataset(std::string const &name, int rows, int cols)
Create a new empty dataset in the file.
void FlushFile()
Flush any data in the HDF5 buffer to the file.
void CreateGroup(std::string const &name)
Create a new group in the file.
bool DoesFileExist(const std::string &name) const
void WriteStringAttribute(std::string const &datasetName, std::string const &attributeName, std::string const &attribute)
Write a string attribute to a dataset or group.
std::string GetStringAttribute(std::string const &datasetName, std::string const &attributeName) const
Read a string attribute from the HDF5 file.
void WriteMatrix(std::string name, Eigen::Matrix< scalarType, fixedRows, fixedCols > const &dataset)
Write data set to file.
bool IsDataSet(std::string const &name) const
void Copy(std::string const &destName, std::shared_ptr< HDF5File > srcFile, std::string const &srcName)
Copy the contents of one dataset into another.
bool IsGroup(std::string const &name) const
std::vector< std::string > GetChildren(std::string base="/") const
Get a list of immediate children of a group.
Eigen::VectorXi GetDataSetSize(std::string const name) const
Get the size of a dataset (rows,cols)
void MergeFile(std::shared_ptr< HDF5File > const &otherFile)
Merge another file into this file.
bool DoesDataSetExist(std::string const &name) const
Check to see if a data set exists.
void WritePartialMatrix(std::string const &name, Eigen::Matrix< scalarType, fixedRows, fixedCols > const &data, unsigned int const row, unsigned int const col)
Write part of a data set to file.
Eigen::Matrix< scalarType, fixedRows, fixedCols > ReadPartialMatrix(std::string const &name, int rowStart, int colStart, int numRows, int numCols) const
void WriteVectorAttribute(std::string const &datasetName, std::string const &attributeName, Eigen::Matrix< scalarType, fixedRows, 1 > const &attribute)
Add a vector attricute to a dataset or group.
Eigen::Matrix< scalarType, fixedRows, fixedCols > ReadMatrix(std::string const &name) const
Read a dataset from the file.
scalarType GetScalarAttribute(std::string const &datasetName, std::string const &attributeName) const
Read a scalar attribute from the HDF5 file.
void WriteScalarAttribute(std::string const &datasetName, std::string const &attributeName, attType const &attribute)
Add a scalar attribute to a dataset or group.
std::string filename
The name of the file.
herr_t H5LTset_attribute(hid_t loc_id, const char *obj_name, const char *attr_name, hid_t mem_type_id, void *buffer, size_t size)
std::string GetParentPath(std::string const &base)