DOLFIN
DOLFIN C++ interface
Public Member Functions | Static Public Member Functions | List of all members
dolfin::HDF5Interface Class Reference

#include <HDF5Interface.h>

Public Member Functions

template<typename T >
void write_dataset (const hid_t file_handle, const std::string dataset_path, const std::vector< T > &data, const std::pair< std::int64_t, std::int64_t > range, const std::vector< int64_t > global_size, bool use_mpi_io, bool use_chunking)
 
template<>
void add_attribute_value (const hid_t dset_id, const std::string attribute_name, const std::string &attribute_value)
 
template<>
void get_attribute_value (const hid_t attr_type, const hid_t attr_id, std::string &attribute_value)
 

Static Public Member Functions

static hid_t open_file (MPI_Comm mpi_comm, const std::string filename, const std::string mode, const bool use_mpi_io)
 Open HDF5 and return file descriptor.
 
static void close_file (const hid_t hdf5_file_handle)
 Close HDF5 file.
 
static void flush_file (const hid_t hdf5_file_handle)
 
static std::string get_filename (hid_t hdf5_file_handle)
 
template<typename T >
static void write_dataset (const hid_t file_handle, const std::string dataset_path, const std::vector< T > &data, const std::pair< std::int64_t, std::int64_t > range, const std::vector< std::int64_t > global_size, bool use_mpio, bool use_chunking)
 
template<typename T >
static void read_dataset (const hid_t file_handle, const std::string dataset_path, const std::pair< std::int64_t, std::int64_t > range, std::vector< T > &data)
 
static bool has_group (const hid_t hdf5_file_handle, const std::string group_name)
 Check for existence of group in HDF5 file.
 
static bool has_dataset (const hid_t hdf5_file_handle, const std::string dataset_path)
 Check for existence of dataset in HDF5 file.
 
static void add_group (const hid_t hdf5_file_handle, const std::string dataset_path)
 Add group to HDF5 file.
 
static int dataset_rank (const hid_t hdf5_file_handle, const std::string dataset_path)
 Get dataset rank.
 
static int num_datasets_in_group (const hid_t hdf5_file_handle, const std::string group_name)
 Return number of data sets in a group.
 
static std::vector< std::int64_t > get_dataset_shape (const hid_t hdf5_file_handle, const std::string dataset_path)
 Get dataset shape (size of each dimension)
 
static std::vector< std::string > dataset_list (const hid_t hdf5_file_handle, const std::string group_name)
 Return list all datasets in named group of file.
 
static const std::string get_attribute_type (const hid_t hdf5_file_handle, const std::string dataset_path, const std::string attribute_name)
 Get type of attribute.
 
template<typename T >
static void get_attribute (const hid_t hdf5_file_handle, const std::string dataset_path, const std::string attribute_name, T &attribute_value)
 Get a named attribute of a dataset of known type.
 
template<typename T >
static void add_attribute (const hid_t hdf5_file_handle, const std::string dataset_path, const std::string attribute_name, const T &attribute_value)
 Add attribute to dataset or group.
 
static void delete_attribute (const hid_t hdf5_file_handle, const std::string dataset_path, const std::string attribute_name)
 Delete an attribute from a dataset or group.
 
static bool has_attribute (const hid_t hdf5_file_handle, const std::string dataset_path, const std::string attribute_name)
 Check if an attribute exists on a dataset or group.
 
static const std::vector< std::string > list_attributes (const hid_t hdf5_file_handle, const std::string dataset_path)
 
static void set_mpi_atomicity (const hid_t hdf5_file_handle, const bool atomic)
 
static bool get_mpi_atomicity (const hid_t hdf5_file_handle)
 

Detailed Description

This class wraps HDF5 function calls. HDF5 function calls should only appear in a member function of this class and not elsewhere in the library.

Member Function Documentation

void HDF5Interface::flush_file ( const hid_t  hdf5_file_handle)
static

Flush data to file to improve data integrity after interruption

bool HDF5Interface::get_mpi_atomicity ( const hid_t  hdf5_file_handle)
static
template<typename T >
void dolfin::HDF5Interface::read_dataset ( const hid_t  file_handle,
const std::string  dataset_path,
const std::pair< std::int64_t, std::int64_t >  range,
std::vector< T > &  data 
)
inlinestatic

Read data from a HDF5 dataset "dataset_path" as defined by range blocks on each process range: the local range on this processor data: a flattened 1D array of values. If range = {-1, -1}, then all data is read on this process.

void HDF5Interface::set_mpi_atomicity ( const hid_t  hdf5_file_handle,
const bool  atomic 
)
static

Set MPI atomicity. See https://support.hdfgroup.org/HDF5/doc/RM/RM_H5F.html#File-SetMpiAtomicity and https://www.open-mpi.org/doc/v2.0/man3/MPI_File_set_atomicity.3.php Writes must be followed by an MPI_Barrier on the communicator before any subsequent reads are guaranteed to return the same data.

template<typename T >
static void dolfin::HDF5Interface::write_dataset ( const hid_t  file_handle,
const std::string  dataset_path,
const std::vector< T > &  data,
const std::pair< std::int64_t, std::int64_t >  range,
const std::vector< std::int64_t >  global_size,
bool  use_mpio,
bool  use_chunking 
)
static

Write data to existing HDF file as defined by range blocks on each process data: data to be written, flattened into 1D vector range: the local range on this processor global_size: the global multidimensional shape of the array use_mpio: whether using MPI or not use_chunking: whether using chunking or not


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