MUQ  0.4.3
muq::Utilities::HDF5File Class Reference

A wrapper around HDF5 (and PHDF5) to handle file I/O. More...

#include <HDF5File.h>

Inheritance diagram for muq::Utilities::HDF5File:

Detailed Description

A wrapper around HDF5 (and PHDF5) to handle file I/O.

Definition at line 23 of file HDF5File.h.

Public Member Functions

 HDF5File (std::string const &filename_)
 Open or create the file. More...
 
virtual ~HDF5File ()
 If HDF5File is destroyed, the file should be closed. More...
 
void Open (std::string const &filename_)
 Opens or creates the file. More...
 
void Close ()
 Close the file. More...
 
void Copy (std::string const &destName, std::shared_ptr< HDF5File > srcFile, std::string const &srcName)
 Copy the contents of one dataset into another. More...
 
template<typename scalarType , int fixedRows, int fixedCols>
void WriteMatrix (std::string name, Eigen::Matrix< scalarType, fixedRows, fixedCols > const &dataset)
 Write data set to file. More...
 
template<typename scalarType , int fixedRows, int fixedCols>
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. More...
 
template<typename scalarType = double, int fixedRows = Eigen::Dynamic, int fixedCols = Eigen::Dynamic>
Eigen::Matrix< scalarType, fixedRows, fixedCols > ReadMatrix (std::string const &name) const
 Read a dataset from the file. More...
 
template<typename scalarType = double, int fixedRows = Eigen::Dynamic, int fixedCols = Eigen::Dynamic>
Eigen::Matrix< scalarType, fixedRows, fixedCols > ReadPartialMatrix (std::string const &name, int rowStart, int colStart, int numRows, int numCols) const
 
bool DoesGroupExist (std::string const &name) const
 Check to see if a group exists. More...
 
bool DoesDataSetExist (std::string const &name) const
 Check to see if a data set exists. More...
 
bool IsDataSet (std::string const &name) const
 
bool IsGroup (std::string const &name) const
 
Eigen::VectorXi GetDataSetSize (std::string const name) const
 Get the size of a dataset (rows,cols) More...
 
void CreateGroup (std::string const &name)
 Create a new group in the file. More...
 
template<typename scalarType >
void CreateDataset (std::string const &name, int rows, int cols)
 Create a new empty dataset in the file. More...
 
template<typename scalarType , int fixedRows>
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. More...
 
template<typename attType = double>
void WriteScalarAttribute (std::string const &datasetName, std::string const &attributeName, attType const &attribute)
 Add a scalar attribute to a dataset or group. More...
 
void WriteStringAttribute (std::string const &datasetName, std::string const &attributeName, std::string const &attribute)
 Write a string attribute to a dataset or group. More...
 
template<typename scalarType >
scalarType GetScalarAttribute (std::string const &datasetName, std::string const &attributeName) const
 Read a scalar attribute from the HDF5 file. More...
 
template<typename scalarType >
Eigen::Matrix< scalarType, Eigen::Dynamic, 1 > GetVectorAttribute (std::string const &datasetName, std::string const &attributeName) const
 Read a vector attribute from the HDF5 file. More...
 
std::string GetStringAttribute (std::string const &datasetName, std::string const &attributeName) const
 Read a string attribute from the HDF5 file. More...
 
std::vector< std::string > GetChildren (std::string base="/") const
 Get a list of immediate children of a group. More...
 
void FlushFile ()
 Flush any data in the HDF5 buffer to the file. More...
 
void MergeFile (std::shared_ptr< HDF5File > const &otherFile)
 Merge another file into this file. More...
 

Public Attributes

hid_t fileID = -1
 The HDF5 file ID. More...
 
std::string filename = ""
 The name of the file. More...
 

Constructor & Destructor Documentation

◆ HDF5File()

HDF5File::HDF5File ( std::string const &  filename_)

Open or create the file.

If the input file exists, it is opened; if the does not exist, it is created.

Parameters
[in]filename_The name (including the path) of the file
[in]writeThe write node that does all of the write-to-file (defaults to 0 — only important if MPI is enabled)
[in]multipleFilesTrue: each rank opens filename seperately and write to its own file (each rank should have a unique filename); False (default): rank 0 opens filename and the other send data to it to write to file.

Definition at line 8 of file HDF5File.cpp.

References fileID, and Open().

◆ ~HDF5File()

HDF5File::~HDF5File ( )
virtual

If HDF5File is destroyed, the file should be closed.

Definition at line 17 of file HDF5File.cpp.

References Close(), and fileID.

Member Function Documentation

◆ Close()

void HDF5File::Close ( )

Close the file.

If the file is not open (HDF5File::fileID is <= 0) then this function does nothing, otherwise the file is closed.

Definition at line 58 of file HDF5File.cpp.

References fileID, filename, and FlushFile().

Referenced by Open(), and ~HDF5File().

◆ Copy()

void HDF5File::Copy ( std::string const &  destName,
std::shared_ptr< HDF5File srcFile,
std::string const &  srcName 
)

Copy the contents of one dataset into another.

Copies the dataset at location srcName in the file srcFile into the location destName of this file. This is a deep copy.

Definition at line 76 of file HDF5File.cpp.

References fileID.

◆ CreateDataset()

template<typename scalarType >
void muq::Utilities::HDF5File::CreateDataset ( std::string const &  name,
int  rows,
int  cols 
)
inline

Create a new empty dataset in the file.

Parameters
[in]nameThe name of the dataset to create. If needed, parent groups will also be created.
[in]rowsThe number of rows in the dataset
[in]colsThe number of columns in the dataset.
[in]singleProcTrue if only a single processor (asumed to be the write parameter) is calling this function (defaults to false, does not matter if MPI is not enabled)

Definition at line 419 of file HDF5File.h.

References CreateGroup(), DoesDataSetExist(), fileID, and muq::Utilities::GetParentPath().

◆ CreateGroup()

void HDF5File::CreateGroup ( std::string const &  name)

Create a new group in the file.

Parameters
[in]nameThe name of the group

Definition at line 167 of file HDF5File.cpp.

References DoesGroupExist(), fileID, FlushFile(), and muq::Utilities::GetParentPath().

Referenced by CreateDataset(), WriteMatrix(), WriteStringAttribute(), and WriteVectorAttribute().

◆ DoesDataSetExist()

bool HDF5File::DoesDataSetExist ( std::string const &  name) const

Check to see if a data set exists.

Parameters
[in]nameThe dataset name (including the path)
Returns
True if the dataset exists, false if not

Definition at line 110 of file HDF5File.cpp.

References DoesGroupExist(), fileID, and muq::Utilities::GetParentPath().

Referenced by CreateDataset(), GetDataSetSize(), GetStringAttribute(), GetVectorAttribute(), IsDataSet(), ReadMatrix(), ReadPartialMatrix(), WriteMatrix(), WritePartialMatrix(), WriteStringAttribute(), and WriteVectorAttribute().

◆ DoesFileExist()

bool HDF5File::DoesFileExist ( const std::string &  name) const
private

Definition at line 25 of file HDF5File.cpp.

Referenced by Open().

◆ DoesGroupExist()

bool HDF5File::DoesGroupExist ( std::string const &  name) const

Check to see if a group exists.

Parameters
[in]nameThe group name (including the path)
Returns
True if the group exists, false if not

Definition at line 93 of file HDF5File.cpp.

References fileID, and muq::Utilities::GetParentPath().

Referenced by CreateGroup(), DoesDataSetExist(), GetChildren(), GetStringAttribute(), GetVectorAttribute(), IsGroup(), WriteStringAttribute(), and WriteVectorAttribute().

◆ FlushFile()

void HDF5File::FlushFile ( )

Flush any data in the HDF5 buffer to the file.

Definition at line 226 of file HDF5File.cpp.

References fileID.

Referenced by Close(), and CreateGroup().

◆ GetChildren()

std::vector< std::string > HDF5File::GetChildren ( std::string  base = "/") const

Get a list of immediate children of a group.

Definition at line 319 of file HDF5File.cpp.

References DoesGroupExist(), fileID, IsDataSet(), and nlohmann::detail::dtoa_impl::len.

◆ GetDataSetSize()

Eigen::VectorXi HDF5File::GetDataSetSize ( std::string const  name) const

Get the size of a dataset (rows,cols)

An empty size zero vector is returned if the dataset could not be found.

Parameters
[in]nameThe name of the data set (include path)

Definition at line 122 of file HDF5File.cpp.

References DoesDataSetExist(), and fileID.

◆ GetScalarAttribute()

template<typename scalarType >
scalarType muq::Utilities::HDF5File::GetScalarAttribute ( std::string const &  datasetName,
std::string const &  attributeName 
) const
inline

Read a scalar attribute from the HDF5 file.

Parameters
[in]datasetNameName of the dataset or group that the attribute is attached to (starting with a "/" and including the path). If the path does not exist, an assert will fail.
[in]attributeNameThe name of the attribute to read. If the attribute does not exist, the HDF5 library may throw an error.

Definition at line 552 of file HDF5File.h.

◆ GetStringAttribute()

std::string HDF5File::GetStringAttribute ( std::string const &  datasetName,
std::string const &  attributeName 
) const

Read a string attribute from the HDF5 file.

Parameters
[in]datasetNameName of the dataset or group that the attribute is attached to (starting with a "/" and including the path). If the path does not exist, an assert will fail.
[in]attributeNameThe name of the attribute to read. If the attribute does not exist, the HDF5 library may throw an error.

Definition at line 205 of file HDF5File.cpp.

References DoesDataSetExist(), DoesGroupExist(), and fileID.

◆ GetVectorAttribute()

template<typename scalarType >
Eigen::Matrix<scalarType, Eigen::Dynamic, 1> muq::Utilities::HDF5File::GetVectorAttribute ( std::string const &  datasetName,
std::string const &  attributeName 
) const
inline

Read a vector attribute from the HDF5 file.

Parameters
[in]datasetNameName of the dataset or group that the attribute is attached to (starting with a "/" and including the path). If the path does not exist, an assert will fail.
[in]attributeNameThe name of the attribute to read. If the attribute does not exist, the HDF5 library may throw an error.

Definition at line 569 of file HDF5File.h.

References DoesDataSetExist(), DoesGroupExist(), and fileID.

◆ IsDataSet()

bool HDF5File::IsDataSet ( std::string const &  name) const

Definition at line 285 of file HDF5File.cpp.

References DoesDataSetExist(), and fileID.

Referenced by GetChildren().

◆ IsGroup()

bool HDF5File::IsGroup ( std::string const &  name) const

Definition at line 302 of file HDF5File.cpp.

References DoesGroupExist(), and fileID.

◆ MergeFile()

void HDF5File::MergeFile ( std::shared_ptr< HDF5File > const &  otherFile)

Merge another file into this file.

Both files must already be open.

Parameters
[in]otherFileThe file to merge into this one.

Definition at line 261 of file HDF5File.cpp.

References CopyObjectToGlobalFile(), and fileID.

◆ Open()

void HDF5File::Open ( std::string const &  filename_)

Opens or creates the file.

If the input file exists, it is opened; if the does not exist, it is created.

Parameters
[in]filename_The name (including the path) of the file
Returns
The file id of the opened or created file

Definition at line 30 of file HDF5File.cpp.

References Close(), DoesFileExist(), fileID, and filename.

Referenced by HDF5File().

◆ ReadMatrix()

template<typename scalarType = double, int fixedRows = Eigen::Dynamic, int fixedCols = Eigen::Dynamic>
Eigen::Matrix<scalarType, fixedRows, fixedCols> muq::Utilities::HDF5File::ReadMatrix ( std::string const &  name) const
inline

Read a dataset from the file.

This function will return an empty matrix if the data set does not exist.

Unless the matrix is full of doubles, the scalarType template parameter should be set to the expected type. Also, fixed size matrices can be returned by setting the fixedRows and/or fixedCols template parameters to something other than Eigen::Dynamic.

Parameters
[in]nameThe name of the dataset.

Definition at line 239 of file HDF5File.h.

References DoesDataSetExist(), and fileID.

◆ ReadPartialMatrix()

template<typename scalarType = double, int fixedRows = Eigen::Dynamic, int fixedCols = Eigen::Dynamic>
Eigen::Matrix<scalarType, fixedRows, fixedCols> muq::Utilities::HDF5File::ReadPartialMatrix ( std::string const &  name,
int  rowStart,
int  colStart,
int  numRows,
int  numCols 
) const
inline

Definition at line 311 of file HDF5File.h.

References DoesDataSetExist(), and fileID.

◆ WriteMatrix()

template<typename scalarType , int fixedRows, int fixedCols>
void muq::Utilities::HDF5File::WriteMatrix ( std::string  name,
Eigen::Matrix< scalarType, fixedRows, fixedCols > const &  dataset 
)
inline

Write data set to file.

Write a new data set to file. If the dataset (and corresponding group) does not exist it will be created. If the dataset does exist it is deleted and overwritten.

In serial (MPI not enabled):

The dataset is written to the user-provided data set name.

In parallel (MPI enabled):

In parallel, multiple datasets are created and (potentially different but the same size) data is written to each one. The names must be unique

Parameters
[in]nameThe name of the data set
[in]datasetThe matrix or vector to save to file (the dataset)
[in]singleProcTrue if only a single processor (asumed to be the write parameter) is calling this function (defaults to false, does not matter if MPI is not enabled)

Definition at line 70 of file HDF5File.h.

References CreateGroup(), DoesDataSetExist(), fileID, and muq::Utilities::GetParentPath().

◆ WritePartialMatrix()

template<typename scalarType , int fixedRows, int fixedCols>
void muq::Utilities::HDF5File::WritePartialMatrix ( std::string const &  name,
Eigen::Matrix< scalarType, fixedRows, fixedCols > const &  data,
unsigned int const  row,
unsigned int const  col 
)
inline

Write part of a data set to file.

Write part of a data set to file. The data set must already exist and be large enough to hold the data given to this function.

Parameters
[in]nameThe name of the data set
[in]dataThe matrix or vector to save to file (part of the dataset)
[in]rowThe row to start writing the data
[in]colThe column to start writting the data
[in]singleProcTrue if only a single processor (asumed to be the write parameter) is calling this function (defaults to false, does not matter if MPI is not enabled)

Definition at line 171 of file HDF5File.h.

References DoesDataSetExist(), and fileID.

◆ WriteScalarAttribute()

template<typename attType = double>
void muq::Utilities::HDF5File::WriteScalarAttribute ( std::string const &  datasetName,
std::string const &  attributeName,
attType const &  attribute 
)
inline

Add a scalar attribute to a dataset or group.

This function adds a scalar valued numeric attribute (such as metadata) to an existing dataset of group. If the given dataset or group does not exist this function creates it (as a group, not a data set).

Parameters
[in]datasetNameThe name of the dataset or group (starting with a "/" and including the path) that you want to attach this attribute to.
[in]attributeNameThe name of the attribute.
[in]attributeThe value of the attribute.
[in]singleProcTrue if only a single processor (asumed to be the write parameter) is calling this function (defaults to false, does not matter if MPI is not enabled)

Definition at line 522 of file HDF5File.h.

◆ WriteStringAttribute()

void HDF5File::WriteStringAttribute ( std::string const &  datasetName,
std::string const &  attributeName,
std::string const &  attribute 
)

Write a string attribute to a dataset or group.

This function adds a scalar valued numeric attribute (such as metadata) to an existing dataset of group. If the given dataset or group does not exist this function creates it (as a group, not a data set).

Parameters
[in]datasetNameThe name of the dataset or group (starting with a "/" and including the path) that you want to attach this attribute to.
[in]attributeNameThe name of the attribute.
[in]attributeThe value of the attribute.
[in]singleProcTrue if only a single processor (asumed to be the write parameter) is calling this function (defaults to false, does not matter if MPI is not enabled)

Definition at line 190 of file HDF5File.cpp.

References CreateGroup(), DoesDataSetExist(), DoesGroupExist(), and fileID.

◆ WriteVectorAttribute()

template<typename scalarType , int fixedRows>
void muq::Utilities::HDF5File::WriteVectorAttribute ( std::string const &  datasetName,
std::string const &  attributeName,
Eigen::Matrix< scalarType, fixedRows, 1 > const &  attribute 
)
inline

Add a vector attricute to a dataset or group.

This function adds a vector valued numeric attribute (such as metadata) to an existing dataset of group. If the given dataset or group does not exist this function creates it (as a group, not a data set).

Parameters
[in]datasetNameThe name of the dataset or group (starting with a "/" and including the path) to attach the attribute to.
[in]attributeNameThe name of the attribute.
[in]attributeThe value of the attribute.
[in]singleProcTrue if only a single processor (asumed to be the write parameter) is calling this function (defaults to false, does not matter if MPI is not enabled)

Definition at line 491 of file HDF5File.h.

References CreateGroup(), DoesDataSetExist(), DoesGroupExist(), fileID, and muq::Utilities::H5LTset_attribute().

Member Data Documentation

◆ fileID

◆ filename

std::string muq::Utilities::HDF5File::filename = ""

The name of the file.

Defaults to an empty string

Definition at line 634 of file HDF5File.h.

Referenced by Close(), and Open().


The documentation for this class was generated from the following files: