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