Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Member Functions | Protected Member Functions | Private Attributes | List of all members
HDF5::DataSet Class Reference

#include <deal.II/base/hdf5.h>

Inheritance diagram for HDF5::DataSet:
[legend]

Public Member Functions

template<typename Container >
Container read ()
 
template<typename Container >
Container read_selection (const std::vector< hsize_t > &coordinates)
 
template<typename Container >
Container read_hyperslab (const std::vector< hsize_t > &offset, const std::vector< hsize_t > &count)
 
template<typename Container >
Container read_hyperslab (const std::vector< hsize_t > &data_dimensions, const std::vector< hsize_t > &offset, const std::vector< hsize_t > &stride, const std::vector< hsize_t > &count, const std::vector< hsize_t > &block)
 
template<typename number >
void read_none ()
 
template<typename Container >
void write (const Container &data)
 
template<typename Container >
void write_selection (const Container &data, const std::vector< hsize_t > &coordinates)
 
template<typename Container >
void write_hyperslab (const Container &data, const std::vector< hsize_t > &offset, const std::vector< hsize_t > &count)
 
template<typename Container >
void write_hyperslab (const Container &data, const std::vector< hsize_t > &data_dimensions, const std::vector< hsize_t > &offset, const std::vector< hsize_t > &stride, const std::vector< hsize_t > &count, const std::vector< hsize_t > &block)
 
template<typename number >
void write_none ()
 
bool get_query_io_mode () const
 
void set_query_io_mode (const bool new_query_io_mode)
 
std::string get_io_mode ()
 
H5D_mpio_actual_io_mode_t get_io_mode_as_hdf5_type ()
 
std::string get_local_no_collective_cause ()
 
uint32_t get_local_no_collective_cause_as_hdf5_type ()
 
std::string get_global_no_collective_cause ()
 
uint32_t get_global_no_collective_cause_as_hdf5_type ()
 
std::vector< hsize_t > get_dimensions () const
 
unsigned int get_size () const
 
unsigned int get_rank () const
 
- Public Member Functions inherited from HDF5::HDF5Object
template<typename T >
get_attribute (const std::string &attr_name) const
 
template<typename T >
void set_attribute (const std::string &attr_name, const T value)
 
std::string get_name () const
 

Protected Member Functions

 DataSet (const std::string &name, const hid_t &parent_group_id, bool mpi)
 
 DataSet (const std::string &name, const hid_t &parent_group_id, const std::vector< hsize_t > &dimensions, const std::shared_ptr< hid_t > &t_type, const bool mpi)
 
- Protected Member Functions inherited from HDF5::HDF5Object
 HDF5Object (const std::string &name, const bool mpi)
 

Private Attributes

unsigned int rank
 
std::vector< hsize_t > dimensions
 
std::shared_ptr< hid_t > dataspace
 
unsigned int size
 
bool query_io_mode
 
H5D_mpio_actual_io_mode_t io_mode
 
uint32_t local_no_collective_cause
 
uint32_t global_no_collective_cause
 

Additional Inherited Members

- Protected Attributes inherited from HDF5::HDF5Object
const std::string name
 
std::shared_ptr< hid_t > hdf5_reference
 
const bool mpi
 

Detailed Description

This class implements an HDF5 DataSet.

Author
Daniel Garcia-Sanchez, 2018, 2019.

Definition at line 392 of file hdf5.h.

Constructor & Destructor Documentation

◆ DataSet() [1/2]

HDF5::DataSet::DataSet ( const std::string &  name,
const hid_t &  parent_group_id,
bool  mpi 
)
protected

Open dataset. This is an internal constructor. The function Group::open_dataset() should be used to open a dataset.

Definition at line 606 of file hdf5.cc.

◆ DataSet() [2/2]

HDF5::DataSet::DataSet ( const std::string &  name,
const hid_t &  parent_group_id,
const std::vector< hsize_t > &  dimensions,
const std::shared_ptr< hid_t > &  t_type,
const bool  mpi 
)
protected

Create dataset. This is an internal constructor. The function Group::create_dataset() should be used to create a dataset.

Definition at line 657 of file hdf5.cc.

Member Function Documentation

◆ read()

template<typename Container >
Container HDF5::DataSet::read ( )

Reads all the data of the dataset.

Datatype conversion takes place at the time of the read operation and is automatic. See the Data Transfer: Datatype Conversion and Selection section in the HDF5 User's Guide.

Container can be std::vector<float>, std::vector<double>, std::vector<std::complex<float>>, std::vector<std::complex<double>>, std::vector<int>, std::vector<unsigned int>, Vector<float>, Vector<double>, Vector<std::complex<float>>, Vector<std::complex<double>>, FullMatrix<float>, FullMatrix<double>, FullMatrix<std::complex<float>> or FullMatrix<std::complex<double>>.

Definition at line 708 of file hdf5.cc.

◆ read_selection()

template<typename Container >
Container HDF5::DataSet::read_selection ( const std::vector< hsize_t > &  coordinates)

Reads data of a subset of the dataset.

Datatype conversion takes place at the time of the read operation and is automatic. See the Data Transfer: Datatype Conversion and Selection section in the HDF5 User's Guide.

The selected elements can be scattered and take any shape in the dataset. For example, in the case of a dataset with rank 4 a selection of 3 points will be described by a 3-by-4 array. Note the indexing is zero-based. To select the points (1,1,1,1), (14,6,12,18), and (8,22,30,22), the point selection array would be as follows:

0 0 0 0
13 5 11 17
7 21 29 21

Parallel HDF5 supports collective I/O on point selections.

Datatype conversion takes place at the time of the read operation and is automatic. See the Data Transfer: Datatype Conversion and Selection section in the HDF5 User's Guide.

Definition at line 743 of file hdf5.cc.

◆ read_hyperslab() [1/2]

template<typename Container >
Container HDF5::DataSet::read_hyperslab ( const std::vector< hsize_t > &  offset,
const std::vector< hsize_t > &  count 
)

Reads a hyperslab from the dataset. The parameters are summarized below:

  • offset: The starting location for the hyperslab.
  • count: The number of elements to select along each dimension.

When reading a hyperslab, HDF5 also allows to provide "stride" and "block" parameters (see the HDF5 documentation). These are not used by the current function and set to nullptr. However these parameters can be used with the function read_hyperslab(const std::vector<hsize_t> &, const std::vector<hsize_t> &, const std::vector<hsize_t> &, const std::vector<hsize_t> &, const std::vector<hsize_t> &)

See the Dataspaces and Data Transfer section in the HDF5 User's Guide. See as well the H5Sselect_hyperslab definition.

Datatype conversion takes place at the time of a read or write and is automatic. See the Data Transfer: Datatype Conversion and Selection section in the HDF5 User's Guide.

Container can be std::vector<float>, std::vector<double>, std::vector<std::complex<float>>, std::vector<std::complex<double>>, std::vector<int>, std::vector<unsigned int>, Vector<float>, Vector<double>, Vector<std::complex<float>>, Vector<std::complex<double>>, FullMatrix<float>, FullMatrix<double>, FullMatrix<std::complex<float>> or FullMatrix<std::complex<double>>.

Definition at line 795 of file hdf5.cc.

◆ read_hyperslab() [2/2]

template<typename Container >
Container HDF5::DataSet::read_hyperslab ( const std::vector< hsize_t > &  data_dimensions,
const std::vector< hsize_t > &  offset,
const std::vector< hsize_t > &  stride,
const std::vector< hsize_t > &  count,
const std::vector< hsize_t > &  block 
)

Writes a data hyperslab to the dataset. The parameters are summarized below:

  • dataset_dimensions: the dimensions of the data memory block.
  • offset: The starting location for the hyperslab.
  • stride: The number of elements to separate each element or block to be selected.
  • count: The number of elements or blocks to select along each dimension.
  • block: The size of the block selected from the dataspace.

See the Dataspaces and Data Transfer section in the HDF5 User's Guide. See as well the H5Sselect_hyperslab definition.

Datatype conversion takes place at the time of a read or write and is automatic. See the Data Transfer: Datatype Conversion and Selection section in the HDF5 User's Guide.

Container can be std::vector<float>, std::vector<double>, std::vector<std::complex<float>>, std::vector<std::complex<double>>, std::vector<int>, std::vector<unsigned int>, Vector<float>, Vector<double>, Vector<std::complex<float>>, Vector<std::complex<double>>, FullMatrix<float>, FullMatrix<double>, FullMatrix<std::complex<float>> or FullMatrix<std::complex<double>>.

Definition at line 849 of file hdf5.cc.

◆ read_none()

template<typename number >
void HDF5::DataSet::read_none ( )

This function does not read any data, but it can contribute to a collective read call. number can be float, double, std::complex<float>, std::complex<double>, int or unsigned int.

Datatype conversion takes place at the time of a read or write and is automatic. See the Data Transfer: Datatype Conversion and Selection section in the HDF5 User's Guide.

Definition at line 902 of file hdf5.cc.

◆ write()

template<typename Container >
void HDF5::DataSet::write ( const Container &  data)

Writes data in the dataset. number can be float, double, std::complex<float>, std::complex<double>, int or unsigned int.

Datatype conversion takes place at the time of a read or write and is automatic. See the Data Transfer: Datatype Conversion and Selection section in the HDF5 User's Guide.

Container can be std::vector<float>, std::vector<double>, std::vector<std::complex<float>>, std::vector<std::complex<double>>, std::vector<int>, std::vector<unsigned int>, Vector<float>, Vector<double>, Vector<std::complex<float>>, Vector<std::complex<double>>, FullMatrix<float>, FullMatrix<double>, FullMatrix<std::complex<float>> or FullMatrix<std::complex<double>>.

Definition at line 942 of file hdf5.cc.

◆ write_selection()

template<typename Container >
void HDF5::DataSet::write_selection ( const Container &  data,
const std::vector< hsize_t > &  coordinates 
)

Writes data to a subset of the dataset. number can be float, double, std::complex<float>, std::complex<double>, int or unsigned int.

The selected elements can be scattered and take any shape in the dataset. For example, in the case of a dataset with rank 4 a selection of 3 points will be described by a 3-by-4 array. Note the indexing is zero-based. To select the points (1,1,1,1), (14,6,12,18), and (8,22,30,22), the point selection array would be as follows:

0 0 0 0
13 5 11 17
7 21 29 21

Parallel HDF5 supports collective I/O on point selections.

Datatype conversion takes place at the time of a read or write and is automatic. See the Data Transfer: Datatype Conversion and Selection section in the HDF5 User's Guide.

Definition at line 974 of file hdf5.cc.

◆ write_hyperslab() [1/2]

template<typename Container >
void HDF5::DataSet::write_hyperslab ( const Container &  data,
const std::vector< hsize_t > &  offset,
const std::vector< hsize_t > &  count 
)

Writes a data hyperslab to the dataset. The parameters are summarized below:

  • offset: The starting location for the hyperslab.
  • count: The number of elements to select along each dimension.

When writing a hyperslab, HDF5 also allows to provide "stride" and "block" parameters (see the HDF5 documentation). These are not used by the current function and set to nullptr. However these parameters can be used with the function write_hyperslab(const Container &data, const std::vector<hsize_t> &data_dimensions, const std::vector<hsize_t> &offset, const std::vector<hsize_t> &stride, const std::vector<hsize_t> &count, const std::vector<hsize_t> &block).

See the Dataspaces and Data Transfer section in the HDF5 User's Guide. See as well the H5Sselect_hyperslab definition.

Datatype conversion takes place at the time of a read or write and is automatic. See the Data Transfer: Datatype Conversion and Selection section in the HDF5 User's Guide.

Definition at line 1023 of file hdf5.cc.

◆ write_hyperslab() [2/2]

template<typename Container >
void HDF5::DataSet::write_hyperslab ( const Container &  data,
const std::vector< hsize_t > &  data_dimensions,
const std::vector< hsize_t > &  offset,
const std::vector< hsize_t > &  stride,
const std::vector< hsize_t > &  count,
const std::vector< hsize_t > &  block 
)

Writes a data hyperslab to the dataset. The parameters are summarized below:

  • dataset_dimensions: the dimensions of the data memory block.
  • offset: The starting location for the hyperslab.
  • stride: The number of elements to separate each element or block to be selected.
  • count: The number of elements or blocks to select along each dimension.
  • block: The size of the block selected from the dataspace.

See the Dataspaces and Data Transfer section in the HDF5 User's Guide. See as well the H5Sselect_hyperslab definition.

Datatype conversion takes place at the time of a read or write and is automatic. See the Data Transfer: Datatype Conversion and Selection section in the HDF5 User's Guide.

Container can be std::vector<float>, std::vector<double>, std::vector<std::complex<float>>, std::vector<std::complex<double>>, std::vector<int>, std::vector<unsigned int>, Vector<float>, Vector<double>, Vector<std::complex<float>>, Vector<std::complex<double>>, FullMatrix<float>, FullMatrix<double>, FullMatrix<std::complex<float>> or FullMatrix<std::complex<double>>.

Definition at line 1080 of file hdf5.cc.

◆ write_none()

template<typename number >
void HDF5::DataSet::write_none ( )

This function does not write any data, but it can contribute to a collective write call. number can be float, double, std::complex<float>, std::complex<double>, int or unsigned int.

Datatype conversion takes place at the time of a read or write and is automatic. See the Data Transfer: Datatype Conversion and Selection section in the HDF5 User's Guide.

Definition at line 1132 of file hdf5.cc.

◆ get_query_io_mode()

bool HDF5::DataSet::get_query_io_mode ( ) const

This function returns the boolean query_io_mode.

In cases where maximum performance has to be achieved, it is important to make sure that all MPI read/write operations are collective. The HDF5 library provides API routines that can be used after the read/write I/O operations to query the I/O mode. If query_io_mode is set to true, then after every read/write operation the deal.II's HDF5 interface calls the routines H5Pget_mpio_actual_io_mode() and H5Pget_mpio_no_collective_cause(). The results are stored in io_mode, local_no_collective_cause and global_no_collective_cause. We suggest to query the I/O mode only in Debug mode because it requires calling additional HDF5 routines.

Definition at line 1279 of file hdf5.cc.

◆ set_query_io_mode()

void HDF5::DataSet::set_query_io_mode ( const bool  new_query_io_mode)

This function sets the boolean query_io_mode.

Definition at line 1171 of file hdf5.cc.

◆ get_io_mode()

std::string HDF5::DataSet::get_io_mode ( )

This function returns the I/O mode that was used on the last parallel I/O call. See H5Pget_mpio_actual_io_mode.

The return value is a std::string and can be

Value Meaning
H5D_MPIO_NO_COLLECTIVE No collective I/O was performed. Collective I/O was not requested or collective I/O isn't possible on this dataset.
H5D_MPIO_CHUNK_INDEPENDENT HDF5 performed chunk collective optimization schemes and each chunk was accessed independently.
H5D_MPIO_CHUNK_COLLECTIVE HDF5 performed chunk collective optimization and each chunk was accessed collectively.
H5D_MPIO_CHUNK_MIXED HDF5 performed chunk collective optimization and some chunks were accessed independently, some collectively.
H5D_MPIO_CONTIGUOUS_COLLECTIVE Collective I/O was performed on a contiguous dataset.

Definition at line 1187 of file hdf5.cc.

◆ get_io_mode_as_hdf5_type()

H5D_mpio_actual_io_mode_t HDF5::DataSet::get_io_mode_as_hdf5_type ( )

This function returns the I/O mode that was used on the last parallel I/O call. See H5Pget_mpio_actual_io_mode. The return type is H5D_mpio_actual_io_mode_t which corresponds to the value returned by H5Pget_mpio_actual_io_mode.

The return value can be

Value Meaning
H5D_MPIO_NO_COLLECTIVE No collective I/O was performed. Collective I/O was not requested or collective I/O isn't possible on this dataset.
H5D_MPIO_CHUNK_INDEPENDENT HDF5 performed chunk collective optimization and each chunk was accessed independently.
H5D_MPIO_CHUNK_COLLECTIVE HDF5 performed chunk collective optimization and each chunk was accessed collectively.
H5D_MPIO_CHUNK_MIXED HDF5 performed chunk collective optimization and some chunks were accessed independently, some collectively.
H5D_MPIO_CONTIGUOUS_COLLECTIVE Collective I/O was performed on a contiguous dataset.

Definition at line 1221 of file hdf5.cc.

◆ get_local_no_collective_cause()

std::string HDF5::DataSet::get_local_no_collective_cause ( )

This function returns the local causes that broke collective I/O on the last parallel I/O call. See H5Pget_mpio_no_collective_cause.

The return value is a string and can be

Value Meaning
H5D_MPIO_COLLECTIVE Collective I/O was performed successfully.
H5D_MPIO_SET_INDEPENDENT Collective I/O was not performed because independent I/O was requested.
H5D_MPIO_DATATYPE_CONVERSION Collective I/O was not performed because datatype conversions were required.
H5D_MPIO_DATA_TRANSFORMS Collective I/O was not performed because data transforms needed to be applied.
H5D_MPIO_SET_MPIPOSIX Collective I/O was not performed because the selected file driver was MPI-POSIX.
H5D_MPIO_NOT_SIMPLE_OR_SCALAR_DATASPACES Collective I/O was not performed because one of the dataspaces was neither simple nor scalar.
H5D_MPIO_POINT_SELECTIONS Collective I/O was not performed because there were point selections in one of the dataspaces.
H5D_MPIO_NOT_CONTIGUOUS_OR_CHUNKED_DATASET Collective I/O was not performed because the dataset was neither contiguous nor chunked.
H5D_MPIO_FILTERS Collective I/O was not performed because filters needed to be applied.

Definition at line 1232 of file hdf5.cc.

◆ get_local_no_collective_cause_as_hdf5_type()

uint32_t HDF5::DataSet::get_local_no_collective_cause_as_hdf5_type ( )

This function returns the local causes that broke collective I/O on the last parallel I/O call. See H5Pget_mpio_no_collective_cause. The return type is uint32_t and corresponds to the value returned by H5Pget_mpio_no_collective_cause.

The return value can be

Value Meaning
H5D_MPIO_COLLECTIVE Collective I/O was performed successfully.
H5D_MPIO_SET_INDEPENDENT Collective I/O was not performed because independent I/O was requested.
H5D_MPIO_DATATYPE_CONVERSION Collective I/O was not performed because datatype conversions were required.
H5D_MPIO_DATA_TRANSFORMS Collective I/O was not performed because data transforms needed to be applied.
H5D_MPIO_SET_MPIPOSIX Collective I/O was not performed because the selected file driver was MPI-POSIX.
H5D_MPIO_NOT_SIMPLE_OR_SCALAR_DATASPACES Collective I/O was not performed because one of the dataspaces was neither simple nor scalar.
H5D_MPIO_POINT_SELECTIONS Collective I/O was not performed because there were point selections in one of the dataspaces.
H5D_MPIO_NOT_CONTIGUOUS_OR_CHUNKED_DATASET Collective I/O was not performed because the dataset was neither contiguous nor chunked.
H5D_MPIO_FILTERS Collective I/O was not performed because filters needed to be applied.

Definition at line 1244 of file hdf5.cc.

◆ get_global_no_collective_cause()

std::string HDF5::DataSet::get_global_no_collective_cause ( )

This function retrieves the global causes that broke collective I/O on the last parallel I/O call. See H5Pget_mpio_no_collective_cause.

The return value is a std::string and can be

Value Meaning
H5D_MPIO_COLLECTIVE Collective I/O was performed successfully.
H5D_MPIO_SET_INDEPENDENT Collective I/O was not performed because independent I/O was requested.
H5D_MPIO_DATATYPE_CONVERSION Collective I/O was not performed because datatype conversions were required.
H5D_MPIO_DATA_TRANSFORMS Collective I/O was not performed because data transforms needed to be applied.
H5D_MPIO_SET_MPIPOSIX Collective I/O was not performed because the selected file driver was MPI-POSIX.
H5D_MPIO_NOT_SIMPLE_OR_SCALAR_DATASPACES Collective I/O was not performed because one of the dataspaces was neither simple nor scalar.
H5D_MPIO_POINT_SELECTIONS Collective I/O was not performed because there were point selections in one of the dataspaces.
H5D_MPIO_NOT_CONTIGUOUS_OR_CHUNKED_DATASET Collective I/O was not performed because the dataset was neither contiguous nor chunked.
H5D_MPIO_FILTERS Collective I/O was not performed because filters needed to be applied.

Definition at line 1256 of file hdf5.cc.

◆ get_global_no_collective_cause_as_hdf5_type()

uint32_t HDF5::DataSet::get_global_no_collective_cause_as_hdf5_type ( )

This function returns the global causes that broke collective I/O on the last parallel I/O call. See H5Pget_mpio_no_collective_cause. The return type is uint32_t and corresponds to the value returned by H5Pget_mpio_no_collective_cause.

The return value can be

Value Meaning
H5D_MPIO_COLLECTIVE Collective I/O was performed successfully.
H5D_MPIO_SET_INDEPENDENT Collective I/O was not performed because independent I/O was requested.
H5D_MPIO_DATATYPE_CONVERSION Collective I/O was not performed because datatype conversions were required.
H5D_MPIO_DATA_TRANSFORMS Collective I/O was not performed because data transforms needed to be applied.
H5D_MPIO_SET_MPIPOSIX Collective I/O was not performed because the selected file driver was MPI-POSIX.
H5D_MPIO_NOT_SIMPLE_OR_SCALAR_DATASPACES Collective I/O was not performed because one of the dataspaces was neither simple nor scalar.
H5D_MPIO_POINT_SELECTIONS Collective I/O was not performed because there were point selections in one of the dataspaces.
H5D_MPIO_NOT_CONTIGUOUS_OR_CHUNKED_DATASET Collective I/O was not performed because the dataset was neither contiguous nor chunked.
H5D_MPIO_FILTERS Collective I/O was not performed because filters needed to be applied.

Definition at line 1268 of file hdf5.cc.

◆ get_dimensions()

std::vector< hsize_t > HDF5::DataSet::get_dimensions ( ) const

This function returns the dimensions of the dataset. The vector dimensions is a one-dimensional array of size rank specifying the size of each dimension of the dataset.

Definition at line 1179 of file hdf5.cc.

◆ get_size()

unsigned int HDF5::DataSet::get_size ( ) const

This function returns the total number of elements in the dataset.

Definition at line 1287 of file hdf5.cc.

◆ get_rank()

unsigned int HDF5::DataSet::get_rank ( ) const

This function returns the rank of the dataset.

Definition at line 1295 of file hdf5.cc.

Member Data Documentation

◆ rank

unsigned int HDF5::DataSet::rank
private

Rank of the DataSet

Definition at line 884 of file hdf5.h.

◆ dimensions

std::vector<hsize_t> HDF5::DataSet::dimensions
private

The vector dimensions is a one-dimensional array of size rank specifying the size of each dimension of the dataset.

Definition at line 890 of file hdf5.h.

◆ dataspace

std::shared_ptr<hid_t> HDF5::DataSet::dataspace
private

HDF5 dataspace identifier.

Definition at line 895 of file hdf5.h.

◆ size

unsigned int HDF5::DataSet::size
private

Total number of elements in the dataset.

Definition at line 900 of file hdf5.h.

◆ query_io_mode

bool HDF5::DataSet::query_io_mode
private

If query_io_mode is set to true, then after every read/write operation the deal.II's HDF5 interface calls the routines H5Pget_mpio_actual_io_mode() and H5Pget_mpio_no_collective_cause(). The results are stored in io_mode, local_no_collective_cause and global_no_collective_cause.

Definition at line 911 of file hdf5.h.

◆ io_mode

H5D_mpio_actual_io_mode_t HDF5::DataSet::io_mode
private

I/O mode that was performed on the last parallel I/O call.

Definition at line 916 of file hdf5.h.

◆ local_no_collective_cause

uint32_t HDF5::DataSet::local_no_collective_cause
private

Local causes that broke collective I/O on the last parallel I/O call. See H5Pget_mpio_no_collective_cause.

Definition at line 923 of file hdf5.h.

◆ global_no_collective_cause

uint32_t HDF5::DataSet::global_no_collective_cause
private

Global causes that broke collective I/O on the last parallel I/O call. See H5Pget_mpio_no_collective_cause.

Definition at line 930 of file hdf5.h.


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