HighFive 3.0.0
HighFive - Header-only C++ HDF5 interface
Loading...
Searching...
No Matches
xtensor.hpp
Go to the documentation of this file.
1#pragma once
2
4#include "H5Exception.hpp"
5
6#include <xtensor/xtensor.hpp>
7#include <xtensor/xarray.hpp>
8#include <xtensor/xadapt.hpp>
9
10namespace HighFive {
11namespace details {
12
13template <class XTensor>
14struct xtensor_get_rank;
15
16template <typename T, size_t N, xt::layout_type L>
17struct xtensor_get_rank<xt::xtensor<T, N, L>> {
18 static constexpr size_t value = N;
19};
20
21template <class EC, size_t N, xt::layout_type L, class Tag>
22struct xtensor_get_rank<xt::xtensor_adaptor<EC, N, L, Tag>> {
23 static constexpr size_t value = N;
24};
25
26template <class Derived, class XTensorType, xt::layout_type L>
27struct xtensor_inspector_base {
28 using type = XTensorType;
29 using value_type = typename type::value_type;
30 using base_type = typename inspector<value_type>::base_type;
31 using hdf5_type = base_type;
32
33 static_assert(std::is_same<value_type, base_type>::value,
34 "HighFive's XTensor support only works for scalar elements.");
35
36 static constexpr bool IsConstExprRowMajor = L == xt::layout_type::row_major;
37 static constexpr bool is_trivially_copyable = IsConstExprRowMajor &&
38 std::is_trivially_copyable<value_type>::value &&
39 inspector<value_type>::is_trivially_copyable;
40
41 static constexpr bool is_trivially_nestable = false;
42
43 static size_t getRank(const type& val) {
44 // Non-scalar elements are not supported.
45 return val.shape().size();
46 }
47
48 static const value_type& getAnyElement(const type& val) {
49 return val.unchecked(0);
50 }
51
52 static value_type& getAnyElement(type& val) {
53 return val.unchecked(0);
54 }
55
56 static std::vector<size_t> getDimensions(const type& val) {
57 auto shape = val.shape();
58 return {shape.begin(), shape.end()};
59 }
60
61 static void prepare(type& val, const std::vector<size_t>& dims) {
62 val.resize(Derived::shapeFromDims(dims));
63 }
64
65 static hdf5_type* data(type& val) {
66 if (!is_trivially_copyable) {
67 throw DataSetException("Invalid used of `inspector<XTensor>::data`.");
68 }
69
70 if (val.size() == 0) {
71 return nullptr;
72 }
73
74 return inspector<value_type>::data(getAnyElement(val));
75 }
76
77 static const hdf5_type* data(const type& val) {
78 if (!is_trivially_copyable) {
79 throw DataSetException("Invalid used of `inspector<XTensor>::data`.");
80 }
81
82 if (val.size() == 0) {
83 return nullptr;
84 }
85
86 return inspector<value_type>::data(getAnyElement(val));
87 }
88
89 static void serialize(const type& val, const std::vector<size_t>& dims, hdf5_type* m) {
90 // since we only support scalar types we know all dims belong to us.
91 size_t size = compute_total_size(dims);
92 xt::adapt(m, size, xt::no_ownership(), dims) = val;
93 }
94
95 static void unserialize(const hdf5_type* vec_align,
96 const std::vector<size_t>& dims,
97 type& val) {
98 // since we only support scalar types we know all dims belong to us.
99 size_t size = compute_total_size(dims);
100 val = xt::adapt(vec_align, size, xt::no_ownership(), dims);
101 }
102};
103
104template <class XTensorType, xt::layout_type L>
105struct xtensor_inspector
106 : public xtensor_inspector_base<xtensor_inspector<XTensorType, L>, XTensorType, L> {
107 private:
108 using super = xtensor_inspector_base<xtensor_inspector<XTensorType, L>, XTensorType, L>;
109
110 public:
111 using type = typename super::type;
112 using value_type = typename super::value_type;
113 using base_type = typename super::base_type;
114 using hdf5_type = typename super::hdf5_type;
115
116 static constexpr size_t ndim = xtensor_get_rank<XTensorType>::value;
117 static constexpr size_t min_ndim = ndim + inspector<value_type>::min_ndim;
118 static constexpr size_t max_ndim = ndim + inspector<value_type>::max_ndim;
119
120 static std::array<size_t, ndim> shapeFromDims(const std::vector<size_t>& dims) {
121 std::array<size_t, ndim> shape;
122 std::copy(dims.cbegin(), dims.cend(), shape.begin());
123 return shape;
124 }
125};
126
127template <class XArrayType, xt::layout_type L>
128struct xarray_inspector
129 : public xtensor_inspector_base<xarray_inspector<XArrayType, L>, XArrayType, L> {
130 private:
131 using super = xtensor_inspector_base<xarray_inspector<XArrayType, L>, XArrayType, L>;
132
133 public:
134 using type = typename super::type;
135 using value_type = typename super::value_type;
136 using base_type = typename super::base_type;
137 using hdf5_type = typename super::hdf5_type;
138
139 static constexpr size_t min_ndim = 0 + inspector<value_type>::min_ndim;
140 static constexpr size_t max_ndim = 1024 + inspector<value_type>::max_ndim;
141
142 static const std::vector<size_t>& shapeFromDims(const std::vector<size_t>& dims) {
143 return dims;
144 }
145};
146
147template <typename T, size_t N, xt::layout_type L>
148struct inspector<xt::xtensor<T, N, L>>: public xtensor_inspector<xt::xtensor<T, N, L>, L> {
149 private:
150 using super = xtensor_inspector<xt::xtensor<T, N, L>, L>;
151
152 public:
153 using type = typename super::type;
154 using value_type = typename super::value_type;
155 using base_type = typename super::base_type;
156 using hdf5_type = typename super::hdf5_type;
157};
158
159template <typename T, xt::layout_type L>
160struct inspector<xt::xarray<T, L>>: public xarray_inspector<xt::xarray<T, L>, L> {
161 private:
162 using super = xarray_inspector<xt::xarray<T, L>, L>;
163
164 public:
165 using type = typename super::type;
166 using value_type = typename super::value_type;
167 using base_type = typename super::base_type;
168 using hdf5_type = typename super::hdf5_type;
169};
170
171template <typename CT, class... S>
172struct inspector<xt::xview<CT, S...>>
173 : public xarray_inspector<xt::xview<CT, S...>, xt::layout_type::any> {
174 private:
175 using super = xarray_inspector<xt::xview<CT, S...>, xt::layout_type::any>;
176
177 public:
178 using type = typename super::type;
179 using value_type = typename super::value_type;
180 using base_type = typename super::base_type;
181 using hdf5_type = typename super::hdf5_type;
182};
183
184
185template <class EC, xt::layout_type L, class SC, class Tag>
186struct inspector<xt::xarray_adaptor<EC, L, SC, Tag>>
187 : public xarray_inspector<xt::xarray_adaptor<EC, L, SC, Tag>, xt::layout_type::any> {
188 private:
189 using super = xarray_inspector<xt::xarray_adaptor<EC, L, SC, Tag>, xt::layout_type::any>;
190
191 public:
192 using type = typename super::type;
193 using value_type = typename super::value_type;
194 using base_type = typename super::base_type;
195 using hdf5_type = typename super::hdf5_type;
196};
197
198template <class EC, size_t N, xt::layout_type L, class Tag>
199struct inspector<xt::xtensor_adaptor<EC, N, L, Tag>>
200 : public xtensor_inspector<xt::xtensor_adaptor<EC, N, L, Tag>, xt::layout_type::any> {
201 private:
202 using super = xtensor_inspector<xt::xtensor_adaptor<EC, N, L, Tag>, xt::layout_type::any>;
203
204 public:
205 using type = typename super::type;
206 using value_type = typename super::value_type;
207 using base_type = typename super::base_type;
208 using hdf5_type = typename super::hdf5_type;
209};
210
211} // namespace details
212} // namespace HighFive
Definition assert_compatible_spaces.hpp:15
size_t compute_total_size(const std::vector< size_t > &dims)
Definition compute_total_size.hpp:10