HighFive 3.0.0
HighFive - Header-only C++ HDF5 interface
Loading...
Searching...
No Matches
H5Easy_Eigen.hpp
Go to the documentation of this file.
1/*
2 * Copyright (c), 2017, Adrien Devresse <adrien.devresse@epfl.ch>
3 *
4 * Distributed under the Boost Software License, Version 1.0.
5 * (See accompanying file LICENSE_1_0.txt or copy at
6 * http://www.boost.org/LICENSE_1_0.txt)
7 *
8 */
9#pragma once
10
11#include "../H5Easy.hpp"
12#include "H5Easy_misc.hpp"
13#include "H5Easy_scalar.hpp"
14
15#ifdef H5_USE_EIGEN
16
17#include "../eigen.hpp"
18
19namespace H5Easy {
20
21namespace detail {
22
23template <typename T>
24struct io_impl<T, typename std::enable_if<std::is_base_of<Eigen::DenseBase<T>, T>::value>::type> {
25 using EigenIndex = Eigen::DenseIndex;
26
27 // When creating a dataset for an Eigen object, the shape of the dataset is
28 // 1D for vectors. (legacy reasons)
29 inline static std::vector<size_t> file_shape(const T& data) {
30 if (std::decay<T>::type::RowsAtCompileTime == 1) {
31 return {static_cast<size_t>(data.cols())};
32 }
33 if (std::decay<T>::type::ColsAtCompileTime == 1) {
34 return {static_cast<size_t>(data.rows())};
35 }
36 return inspector<T>::getDimensions(data);
37 }
38
39 // The shape of an Eigen object as used in HighFive core.
40 inline static std::vector<size_t> mem_shape(const T& data) {
41 return inspector<T>::getDimensions(data);
42 }
43
44 // The shape of an Eigen object as used in HighFive core.
45 template <class D>
46 inline static std::vector<size_t> mem_shape(const File& file,
47 const std::string& path,
48 const D& dataset) {
49 std::vector<size_t> dims = dataset.getDimensions();
50
51 if (dims.size() == 1 && T::RowsAtCompileTime == 1) {
52 return std::vector<size_t>{1, dims[0]};
53 }
54 if (dims.size() == 1 && T::ColsAtCompileTime == 1) {
55 return std::vector<size_t>{dims[0], 1};
56 }
57 if (dims.size() == 2) {
58 return dims;
59 }
60
61 throw detail::error(file, path, "H5Easy::load: Inconsistent rank");
62 }
63
64 inline static DataSet dump(File& file,
65 const std::string& path,
66 const T& data,
67 const DumpOptions& options) {
68 using value_type = typename std::decay<T>::type::Scalar;
69
70 std::vector<size_t> file_dims = file_shape(data);
71 std::vector<size_t> mem_dims = mem_shape(data);
72 DataSet dataset = initDataset<value_type>(file, path, file_dims, options);
73 dataset.reshapeMemSpace(mem_dims).write(data);
74 if (options.flush()) {
75 file.flush();
76 }
77 return dataset;
78 }
79
80 inline static T load(const File& file, const std::string& path) {
81 DataSet dataset = file.getDataSet(path);
82 std::vector<size_t> dims = mem_shape(file, path, dataset);
83 return dataset.reshapeMemSpace(dims).template read<T>();
84 }
85
86 inline static Attribute dumpAttribute(File& file,
87 const std::string& path,
88 const std::string& key,
89 const T& data,
90 const DumpOptions& options) {
91 using value_type = typename std::decay<T>::type::Scalar;
92
93 std::vector<size_t> file_dims = file_shape(data);
94 std::vector<size_t> mem_dims = mem_shape(data);
95 Attribute attribute = initAttribute<value_type>(file, path, key, file_dims, options);
96 attribute.reshapeMemSpace(mem_dims).write(data);
97 if (options.flush()) {
98 file.flush();
99 }
100 return attribute;
101 }
102
103 inline static T loadAttribute(const File& file,
104 const std::string& path,
105 const std::string& key) {
106 DataSet dataset = file.getDataSet(path);
107 Attribute attribute = dataset.getAttribute(key);
108 DataSpace dataspace = attribute.getSpace();
109 std::vector<size_t> dims = mem_shape(file, path, dataspace);
110 return attribute.reshapeMemSpace(dims).template read<T>();
111 }
112};
113
114} // namespace detail
115} // namespace H5Easy
116
117#endif // H5_USE_EIGEN
Read/dump DataSets or Attribute using a minimalistic syntax. To this end, the functions are templated...
Definition H5Easy.hpp:62
DataSet dump(File &file, const std::string &path, const T &data, DumpMode mode=DumpMode::Create)
Write object (templated) to a (new) DataSet in an open HDF5 file.
Definition H5Easy_public.hpp:99
T loadAttribute(const File &file, const std::string &path, const std::string &key)
Load a Attribute in an open HDF5 file to an object (templated).
Definition H5Easy_public.hpp:166
Attribute dumpAttribute(File &file, const std::string &path, const std::string &key, const T &data, DumpMode mode=DumpMode::Create)
Write object (templated) to a (new) Attribute in an open HDF5 file.
Definition H5Easy_public.hpp:148
T load(const File &file, const std::string &path, const std::vector< size_t > &idx)
Load entry {i, j, ...} from a DataSet in an open HDF5 file to a scalar.
Definition H5Easy_public.hpp:138