35#ifndef NEKTAR_LIB_UTILITIES_BASIC_UTILS_H5_H
36#define NEKTAR_LIB_UTILITIES_BASIC_UTILS_H5_H
50#define H5_CONSTRUCT(ans, func, args) \
52 hid_t ret = func args; \
54 ErrorUtil::Error(ErrorUtil::efatal, __FILE__, __LINE__, \
55 "HDF5 error in API function " #func, 0); \
58#define H5_CALL(func, args) \
60 herr_t ret = func args; \
62 ErrorUtil::Error(ErrorUtil::efatal, __FILE__, __LINE__, \
63 "HDF5 error in API function " #func, 0); \
66class Error :
public std::exception
95class Object :
public std::enable_shared_from_this<Object>
108 inline operator hid_t()
const
162 void SetChunk(
const std::vector<hsize_t> &dims);
188 return !(*
this == other);
201 static herr_t
helper(hid_t g_id,
const char *
name,
202 const H5L_info_t *info,
void *op_data);
228 const std::string &
name,
const std::vector<T> &data,
240 const std::string &
name,
271 return !(*
this == other);
279 static herr_t
helper(hid_t g_id,
const char *
name,
280 const H5A_info_t *info,
void *op_data);
316 DataSpace(hsize_t size, hsize_t max = H5S_UNLIMITED - 1);
317 DataSpace(
const std::vector<hsize_t> &dims);
318 DataSpace(
const std::vector<hsize_t> &dims,
319 const std::vector<hsize_t> &max_dims);
322 void SelectRange(
const hsize_t start,
const hsize_t count);
323 void AppendRange(
const hsize_t start,
const hsize_t count);
326 const std::vector<hsize_t> count);
328 const std::vector<hsize_t> count);
331 const std::vector<hsize_t> &coords);
336 std::vector<hsize_t>
GetDims();
413 void Add(std::string
name,
size_t offset, hid_t type)
415 H5Tinsert(
m_Id,
name.c_str(), offset, type);
419 hid_t strtype = H5Tcopy(H5T_C_S1);
420 H5Tset_size(strtype, size);
421 H5Tinsert(
m_Id,
name.c_str(), offset, strtype);
509 template <
class T>
void Write(
const std::vector<T> &data)
512 H5_CALL(H5Dwrite, (
m_Id, mem_t->GetId(), H5S_ALL, H5S_ALL, H5P_DEFAULT,
522 H5_CALL(H5Dwrite, (
m_Id, mem_t->GetId(), memspace->GetId(),
523 filespace->GetId(), dxpl->GetId(), &data[0]));
532 H5Dwrite(
m_Id, dt->GetId(), memspace->GetId(), filespace->GetId(),
533 dxpl->GetId(), &data[0]);
540 H5_CALL(H5Dwrite, (
m_Id, type->GetId(), H5S_ALL, filespace->GetId(),
541 dxpl->GetId(), s.c_str()));
549 ret =
new const char *[s.size()];
550 for (
size_t i = 0; i < s.size(); ++i)
552 ret[i] = s[i].c_str();
554 H5Dwrite(
m_Id, type->GetId(), H5S_ALL, filespace->GetId(),
558 template <
class T>
void Read(std::vector<T> &data)
562 ASSERTL0(H5Sget_simple_extent_ndims(space->GetId()) == 1,
563 "vector data not 1D");
565 H5Sget_simple_extent_dims(space->GetId(), &len, &maxdim);
569 H5_CALL(H5Dread, (
m_Id, mem_t->GetId(), H5S_ALL, H5S_ALL, H5P_DEFAULT,
578 len = H5Sget_select_npoints(filespace->GetId());
583 H5_CALL(H5Dread, (
m_Id, mem_t->GetId(), memspace->GetId(),
584 filespace->GetId(), dxpl->GetId(), &data[0]));
588 std::vector<std::vector<int>> &coords,
594 int w = coords[0].size();
596 hsize_t *cds =
new hsize_t[coords.size() *
w];
597 for (
int i = 0; i < coords.size(); i++)
599 for (
int j = 0; j < coords[i].size(); j++)
601 cds[i *
w + j] = hsize_t(coords[i][j]);
605 H5Sselect_elements(filespace->GetId(), H5S_SELECT_SET, coords.size(),
610 len = H5Sget_select_npoints(filespace->GetId());
615 H5_CALL(H5Dread, (
m_Id, mem_t->GetId(), memspace->GetId(),
616 filespace->GetId(), dxpl->GetId(), &data[0]));
618 H5Sselect_all(filespace->GetId());
627 std::vector<hsize_t> dims = filespace->GetDims();
628 rdata = (
char **)malloc(dims[0] *
sizeof(
char *));
630 H5_CALL(H5Dread, (
m_Id, tp->GetId(), H5S_ALL, filespace->GetId(),
631 dxpl->GetId(), rdata));
633 for (
int i = 0; i < dims[0]; i++)
635 data.push_back(std::string(rdata[i]));
649 T>::Convert(
const T &obj)
663 return PredefinedDataType::Native<T>();
683 return Converter::Convert(obj);
690 return Converter::Deconvert(obj);
703 return std::string(obj);
727 H5_CALL(H5Awrite, (attr->GetId(), type->GetId(),
733 const std::vector<T> &value)
743 const void *converted_buf =
nullptr;
746 converted_vals.resize(value.size());
747 for (
size_t i = 0; i < value.size(); ++i)
751 converted_buf = &converted_vals[0];
755 converted_buf = &value[0];
758 H5_CALL(H5Awrite, (attr->GetId(), type->GetId(), converted_buf));
770 H5_CALL(H5Aread, (attr->GetId(), type->GetId(),
777 std::vector<T> &value)
786 ASSERTL0(H5Sget_simple_extent_ndims(space->GetId()) == 1,
787 "vector data not 1D");
789 H5Sget_simple_extent_dims(space->GetId(), &len, &maxdim);
792 void *converted_buf =
nullptr;
795 converted_vals.resize(len);
796 converted_buf = &converted_vals[0];
800 converted_buf = &value[0];
803 H5_CALL(H5Aread, (attr->GetId(), type->GetId(), converted_buf));
807 typename Vec::iterator src = converted_vals.begin(),
808 end = converted_vals.end();
809 typename std::vector<T>::iterator dest = value.begin();
810 for (; src != end; ++src, ++dest)
819 const std::string &
name,
const std::vector<T> &data,
826 dataset->Write(data);
#define ASSERTL0(condition, msg)
#define H5_CALL(func, args)
static AttributeSharedPtr Open(hid_t parent, const std::string &name)
DataSpaceSharedPtr GetSpace() const
static AttributeSharedPtr Create(hid_t parent, const std::string &name, DataTypeSharedPtr type, DataSpaceSharedPtr space)
bool operator==(const AttrIterator &other) const
static herr_t helper(hid_t g_id, const char *name, const H5A_info_t *info, void *op_data)
const std::string & operator*()
AttrIterator(CanHaveAttributesSharedPtr obj, hsize_t idx=0)
bool operator!=(const AttrIterator &other) const
CanHaveAttributesSharedPtr m_obj
AttrIterator & operator++()
std::string m_currentName
Mixin for objects that can have attributes (Group, DataSet, DataType)
AttributeSharedPtr OpenAttribute(const std::string &name)
void GetAttribute(const std::string &name, T &value)
AttrIterator attr_begin()
AttributeSharedPtr CreateAttribute(const std::string &name, DataTypeSharedPtr type, DataSpaceSharedPtr space)
void SetAttribute(const std::string &name, const T &value)
LinkIterator & operator++()
static herr_t helper(hid_t g_id, const char *name, const H5L_info_t *info, void *op_data)
std::string m_currentName
bool operator!=(const LinkIterator &other) const
LinkIterator(CanHaveGroupsDataSetsSharedPtr grp, hsize_t idx=0)
std::string GetName() const
bool operator==(const LinkIterator &other) const
CanHaveGroupsDataSetsSharedPtr m_grp
const std::string & operator*()
Mixin for objects that contain groups and datasets (Group and File)
virtual hsize_t v_GetNumElements()=0
GroupSharedPtr OpenGroup(const std::string &name, PListSharedPtr accessPL=PList::Default()) const
DataSetSharedPtr OpenDataSet(const std::string &name, PListSharedPtr accessPL=PList::Default()) const
friend class key_iterator
GroupSharedPtr CreateGroup(const std::string &name, PListSharedPtr createPL=PList::Default(), PListSharedPtr accessPL=PList::Default())
DataSetSharedPtr CreateDataSet(const std::string &name, DataTypeSharedPtr type, DataSpaceSharedPtr space, PListSharedPtr createPL=PList::Default(), PListSharedPtr accessPL=PList::Default())
DataSetSharedPtr CreateWriteDataSet(const std::string &name, const std::vector< T > &data, PListSharedPtr createPL=PList::Default(), PListSharedPtr accessPL=PList::Default())
bool ContainsDataSet(std::string nm)
void AddString(std::string name, size_t offset, size_t size)
static CompoundDataTypeSharedPtr Create(size_t sz)
void Add(std::string name, size_t offset, hid_t type)
void Read(std::vector< T > &data)
void Write(const std::vector< T > &data, DataSpaceSharedPtr filespace, PListSharedPtr dxpl=PList::Default())
void ReadVectorString(std::vector< std::string > &data, DataSpaceSharedPtr filespace, PListSharedPtr dxpl=PList::Default())
void WriteString(std::string s, DataSpaceSharedPtr filespace, DataTypeSharedPtr type, PListSharedPtr dxpl=PList::Default())
void Read(std::vector< T > &data, DataSpaceSharedPtr filespace, std::vector< std::vector< int > > &coords, PListSharedPtr dxpl=PList::Default())
DataSpaceSharedPtr GetSpace() const
void Read(std::vector< T > &data, DataSpaceSharedPtr filespace, PListSharedPtr dxpl=PList::Default())
void WriteVectorString(std::vector< std::string > s, DataSpaceSharedPtr filespace, DataTypeSharedPtr type, PListSharedPtr dxpl=PList::Default())
void Write(const std::vector< T > &data, DataSpaceSharedPtr filespace, DataTypeSharedPtr dt, PListSharedPtr dxpl=PList::Default())
void Write(const std::vector< T > &data)
void SelectRange(const hsize_t start, const hsize_t count)
void AppendRange(const hsize_t start, const hsize_t count)
static DataSpaceSharedPtr Null()
static DataSpaceSharedPtr Scalar()
std::vector< hsize_t > GetDims()
static DataSpaceSharedPtr OneD(hsize_t size)
void SetSelection(const hsize_t num_elmt, const std::vector< hsize_t > &coords)
Wrap and HDF5 data type object. Technically this can have attributes, but not really bothered.
static DataTypeSharedPtr OfObject(const T &obj)
DataTypeSharedPtr Copy() const
static DataTypeSharedPtr String(size_t len=0)
static FileSharedPtr Open(const std::string &filename, unsigned mode, PListSharedPtr accessPL=PList::Default())
static FileSharedPtr Create(const std::string &filename, unsigned mode, PListSharedPtr createPL=PList::Default(), PListSharedPtr accessPL=PList::Default())
hsize_t v_GetNumElements() override
std::vector< std::string > GetElementNames()
hsize_t v_GetNumElements() override
CanHaveAttributesSharedPtr operator[](hsize_t idx)
CanHaveAttributesSharedPtr operator[](const std::string &key)
static PListSharedPtr ObjectCreate()
Properties for object creation.
static PListSharedPtr DatasetXfer()
Properties for raw data transfer.
static PListSharedPtr StringCreate()
Properties for character encoding when encoding strings or object names.
static PListSharedPtr FileAccess()
Properties for file access.
void SetDxMpioIndependent()
void SetMpio(CommSharedPtr comm)
void SetDxMpioCollective()
void SetDeflate(const unsigned level=1)
static PListSharedPtr FileCreate()
Properties for file creation.
static PListSharedPtr GroupAccess()
Properties for group access.
static PListSharedPtr DatasetAccess()
Properties for dataset access.
static PListSharedPtr Default()
Default options.
static PListSharedPtr FileMount()
Properties for file mounting.
static PListSharedPtr ObjectCopy()
Properties governing the object copying process.
static PListSharedPtr DatasetCreate()
Properties for dataset creation.
static PListSharedPtr GroupCreate()
Properties for group creation.
void SetChunk(const std::vector< hsize_t > &dims)
static PListSharedPtr AttributeCreate()
Properties for attribute creation.
static PListSharedPtr LinkAccess()
Properties governing link traversal when accessing objects.
static PListSharedPtr LinkCreate()
Properties governing link creation.
static PListSharedPtr DatatypeCreate()
Properties for datatype creation.
static PListSharedPtr DatatypeAccess()
Properties for datatype access.
Predefined HDF data types that must not be closed when done with.
PredefinedDataType(hid_t)
static DataTypeSharedPtr Native()
static DataTypeSharedPtr CS1()
std::shared_ptr< DataSpace > DataSpaceSharedPtr
std::shared_ptr< PList > PListSharedPtr
std::shared_ptr< Object > ObjectSharedPtr
std::shared_ptr< File > FileSharedPtr
std::shared_ptr< DataType > DataTypeSharedPtr
std::shared_ptr< CanHaveGroupsDataSets > CanHaveGroupsDataSetsSharedPtr
std::shared_ptr< Attribute > AttributeSharedPtr
std::shared_ptr< CompoundDataType > CompoundDataTypeSharedPtr
std::shared_ptr< CanHaveAttributes > CanHaveAttributesSharedPtr
std::shared_ptr< Group > GroupSharedPtr
std::shared_ptr< DataSet > DataSetSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
std::vector< double > w(NPUPPER)
const char * ConvertedVectorElemType
static std::string Deconvert(const ConvertedType &obj)
const char * ConvertedType
static ConvertedType Convert(const std::string &obj)
static ConvertedType Convert(const T &obj)
T ConvertedVectorElemType
static T Deconvert(const ConvertedType &obj)
static const bool MustConvert
Traits class for HDF5 data types.
DataTypeConversionPolicy< T > Converter
Converter::ConvertedType ConvertedType
static const void * GetAddress(const ConvertedType &obj)
static void * GetAddress(ConvertedType &obj)
static DataTypeSharedPtr GetType()
static ConvertedType Convert(const T &obj)
static const hid_t NativeType
static DataTypeSharedPtr GetType(const T &obj)
static T Deconvert(const ConvertedType &obj)