HighFive 2.9.0
HighFive - Header-only C++ HDF5 interface
Loading...
Searching...
No Matches
H5Slice_traits_misc.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 <algorithm>
12#include <cassert>
13#include <functional>
14#include <numeric>
15#include <sstream>
16#include <string>
17
18#include "h5d_wrapper.hpp"
19#include "h5s_wrapper.hpp"
20
21#include "H5ReadWrite_misc.hpp"
22#include "H5Converter_misc.hpp"
23
24namespace HighFive {
25
26namespace details {
27
28// map the correct reference to the dataset depending of the layout
29// dataset -> itself
30// subselection -> parent dataset
31inline const DataSet& get_dataset(const Selection& sel) {
32 return sel.getDataset();
33}
34
35inline const DataSet& get_dataset(const DataSet& ds) {
36 return ds;
37}
38
39// map the correct memspace identifier depending of the layout
40// dataset -> entire memspace
41// selection -> resolve space id
42inline hid_t get_memspace_id(const Selection& ptr) {
43 return ptr.getMemSpace().getId();
44}
45
46inline hid_t get_memspace_id(const DataSet&) {
47 return H5S_ALL;
48}
49} // namespace details
50
51inline ElementSet::ElementSet(std::initializer_list<std::size_t> list)
52 : _ids(list) {}
53
54inline ElementSet::ElementSet(std::initializer_list<std::vector<std::size_t>> list)
55 : ElementSet(std::vector<std::vector<std::size_t>>(list)) {}
56
57inline ElementSet::ElementSet(const std::vector<std::size_t>& element_ids)
58 : _ids(element_ids) {}
59
60inline ElementSet::ElementSet(const std::vector<std::vector<std::size_t>>& element_ids) {
61 for (const auto& vec: element_ids) {
62 std::copy(vec.begin(), vec.end(), std::back_inserter(_ids));
63 }
64}
65
66template <typename Derivate>
68 const DataSpace& memspace) const {
69 // Note: The current limitation are that memspace must describe a
70 // packed memspace.
71 //
72 // The reason for this is that we're unable to unpack general
73 // hyperslabs when the memory is not contiguous, e.g.
74 // `std::vector<std::vector<double>>`.
75 const auto& slice = static_cast<const Derivate&>(*this);
76 auto filespace = hyperslab.apply(slice.getSpace());
77
78 return detail::make_selection(memspace, filespace, details::get_dataset(slice));
79}
80
81template <typename Derivate>
82inline Selection SliceTraits<Derivate>::select(const HyperSlab& hyper_slab) const {
83 const auto& slice = static_cast<const Derivate&>(*this);
84 auto filespace = slice.getSpace();
85 filespace = hyper_slab.apply(filespace);
86
87 auto n_elements = detail::h5s_get_select_npoints(filespace.getId());
88 auto memspace = DataSpace(std::array<size_t, 1>{size_t(n_elements)});
89
90 return detail::make_selection(memspace, filespace, details::get_dataset(slice));
91}
92
93
94template <typename Derivate>
95inline Selection SliceTraits<Derivate>::select(const std::vector<size_t>& offset,
96 const std::vector<size_t>& count,
97 const std::vector<size_t>& stride,
98 const std::vector<size_t>& block) const {
99 auto slab = HyperSlab(RegularHyperSlab(offset, count, stride, block));
100 auto memspace = DataSpace(count);
101 return select(slab, memspace);
102}
103
104template <typename Derivate>
105inline Selection SliceTraits<Derivate>::select(const std::vector<size_t>& columns) const {
106 const auto& slice = static_cast<const Derivate&>(*this);
107 const DataSpace& space = slice.getSpace();
108 std::vector<size_t> dims = space.getDimensions();
109
110 std::vector<size_t> counts = dims;
111 counts.back() = 1;
112
113 std::vector<size_t> offsets(dims.size(), 0);
114
115 HyperSlab slab;
116 for (const auto& column: columns) {
117 offsets.back() = column;
118 slab |= RegularHyperSlab(offsets, counts);
119 }
120
121 std::vector<size_t> memdims = dims;
122 memdims.back() = columns.size();
123
124 return select(slab, DataSpace(memdims));
125}
126
127template <typename Derivate>
129 const auto& slice = static_cast<const Derivate&>(*this);
130 const hsize_t* data = nullptr;
131 const DataSpace space = slice.getSpace().clone();
132 const std::size_t length = elements._ids.size();
133 if (length % space.getNumberDimensions() != 0) {
134 throw DataSpaceException(
135 "Number of coordinates in elements picking "
136 "should be a multiple of the dimensions.");
137 }
138 const std::size_t num_elements = length / space.getNumberDimensions();
139 std::vector<hsize_t> raw_elements;
140
141 // optimised at compile time
142 // switch for data conversion on 32bits platforms
143 if (std::is_same<std::size_t, hsize_t>::value) {
144 // `if constexpr` can't be used, thus a reinterpret_cast is needed.
145 data = reinterpret_cast<const hsize_t*>(&(elements._ids[0]));
146 } else {
147 raw_elements.resize(length);
148 std::copy(elements._ids.begin(), elements._ids.end(), raw_elements.begin());
149 data = raw_elements.data();
150 }
151
152 detail::h5s_select_elements(space.getId(), H5S_SELECT_SET, num_elements, data);
153
154 return detail::make_selection(DataSpace(num_elements), space, details::get_dataset(slice));
155}
156
157
158template <typename Derivate>
159template <typename T>
160inline T SliceTraits<Derivate>::read(const DataTransferProps& xfer_props) const {
161 T array;
162 read(array, xfer_props);
163 return array;
164}
165
166
167template <typename Derivate>
168template <typename T>
169inline void SliceTraits<Derivate>::read(T& array, const DataTransferProps& xfer_props) const {
170 const auto& slice = static_cast<const Derivate&>(*this);
171 const DataSpace& mem_space = slice.getMemSpace();
172
173 auto file_datatype = slice.getDataType();
174
175 const details::BufferInfo<T> buffer_info(
176 file_datatype,
177 [&slice]() -> std::string { return details::get_dataset(slice).getPath(); },
178 details::BufferInfo<T>::Operation::read);
179
180 if (!details::checkDimensions(mem_space, buffer_info.n_dimensions)) {
181 std::ostringstream ss;
182 ss << "Impossible to read DataSet of dimensions " << mem_space.getNumberDimensions()
183 << " into arrays of dimensions " << buffer_info.n_dimensions;
184 throw DataSpaceException(ss.str());
185 }
186 auto dims = mem_space.getDimensions();
187
188 auto r = details::data_converter::get_reader<T>(dims, array, file_datatype);
189 read_raw(r.getPointer(), buffer_info.data_type, xfer_props);
190 // re-arrange results
191 r.unserialize(array);
192
193 auto t = buffer_info.data_type;
194 auto c = t.getClass();
195 if (c == DataTypeClass::VarLen || t.isVariableStr()) {
196#if H5_VERSION_GE(1, 12, 0)
197 // This one have been created in 1.12.0
198 (void)
199 detail::h5t_reclaim(t.getId(), mem_space.getId(), xfer_props.getId(), r.getPointer());
200#else
201 // This one is deprecated since 1.12.0
202 (void) detail::h5d_vlen_reclaim(t.getId(),
203 mem_space.getId(),
204 xfer_props.getId(),
205 r.getPointer());
206#endif
207 }
208}
209
210template <typename Derivate>
211template <typename T>
212inline void SliceTraits<Derivate>::read(T* array,
213 const DataType& mem_datatype,
214 const DataTransferProps& xfer_props) const {
215 read_raw(array, mem_datatype, xfer_props);
216}
217
218template <typename Derivate>
219template <typename T>
220inline void SliceTraits<Derivate>::read(T* array, const DataTransferProps& xfer_props) const {
221 read_raw(array, xfer_props);
222}
223
224
225template <typename Derivate>
226template <typename T>
228 const DataType& mem_datatype,
229 const DataTransferProps& xfer_props) const {
230 static_assert(!std::is_const<T>::value,
231 "read() requires a non-const structure to read data into");
232
233 const auto& slice = static_cast<const Derivate&>(*this);
234
235 detail::h5d_read(details::get_dataset(slice).getId(),
236 mem_datatype.getId(),
237 details::get_memspace_id(slice),
238 slice.getSpace().getId(),
239 xfer_props.getId(),
240 static_cast<void*>(array));
241}
242
243
244template <typename Derivate>
245template <typename T>
246inline void SliceTraits<Derivate>::read_raw(T* array, const DataTransferProps& xfer_props) const {
247 using element_type = typename details::inspector<T>::base_type;
248 const DataType& mem_datatype = create_and_check_datatype<element_type>();
249
250 read_raw(array, mem_datatype, xfer_props);
251}
252
253
254template <typename Derivate>
255template <typename T>
256inline void SliceTraits<Derivate>::write(const T& buffer, const DataTransferProps& xfer_props) {
257 const auto& slice = static_cast<const Derivate&>(*this);
258 const DataSpace& mem_space = slice.getMemSpace();
259
260 auto file_datatype = slice.getDataType();
261
262 const details::BufferInfo<T> buffer_info(
263 file_datatype,
264 [&slice]() -> std::string { return details::get_dataset(slice).getPath(); },
265 details::BufferInfo<T>::Operation::write);
266
267 if (!details::checkDimensions(mem_space, buffer_info.n_dimensions)) {
268 std::ostringstream ss;
269 ss << "Impossible to write buffer of dimensions "
270 << details::format_vector(mem_space.getDimensions())
271 << " into dataset with n = " << buffer_info.n_dimensions << " dimensions.";
272 throw DataSpaceException(ss.str());
274 auto w = details::data_converter::serialize<T>(buffer, file_datatype);
275 write_raw(w.getPointer(), buffer_info.data_type, xfer_props);
276}
277
278
279template <typename Derivate>
280template <typename T>
281inline void SliceTraits<Derivate>::write_raw(const T* buffer,
282 const DataType& mem_datatype,
283 const DataTransferProps& xfer_props) {
284 const auto& slice = static_cast<const Derivate&>(*this);
285
286 detail::h5d_write(details::get_dataset(slice).getId(),
287 mem_datatype.getId(),
288 details::get_memspace_id(slice),
289 slice.getSpace().getId(),
290 xfer_props.getId(),
291 static_cast<const void*>(buffer));
292}
293
294
295template <typename Derivate>
296template <typename T>
297inline void SliceTraits<Derivate>::write_raw(const T* buffer, const DataTransferProps& xfer_props) {
298 using element_type = typename details::inspector<T>::base_type;
299 const auto& mem_datatype = create_and_check_datatype<element_type>();
300
301 write_raw(buffer, mem_datatype, xfer_props);
302}
303
305} // namespace HighFive
Exception specific to HighFive DataSpace interface.
Definition H5Exception.hpp:112
Class representing the space (dimensions) of a DataSet.
Definition H5DataSpace.hpp:31
size_t getNumberDimensions() const
Returns the number of dimensions of a DataSpace.
Definition H5Dataspace_misc.hpp:94
std::vector< size_t > getDimensions() const
Returns the size of the dataset in each dimension.
Definition H5Dataspace_misc.hpp:98
DataSpace clone() const
Create a copy of the DataSpace which will have different id.
Definition H5Dataspace_misc.hpp:88
HDF5 Data Type.
Definition H5DataType.hpp:61
Definition H5Slice_traits.hpp:22
ElementSet(std::initializer_list< std::size_t > list)
Create a list of points of N-dimension for selection.
Definition H5Slice_traits_misc.hpp:51
Definition H5Slice_traits.hpp:121
DataSpace apply(const DataSpace &space_) const
Definition H5Slice_traits.hpp:174
hid_t getId() const noexcept
getId
Definition H5Object_misc.hpp:69
HDF5 property Lists.
Definition H5PropertyList.hpp:160
Selection: represent a view on a slice/part of a dataset.
Definition H5Selection.hpp:27
DataSpace getSpace() const noexcept
getSpace
Definition H5Selection_misc.hpp:20
T read(const DataTransferProps &xfer_props=DataTransferProps()) const
Definition H5Slice_traits_misc.hpp:160
void write(const T &buffer, const DataTransferProps &xfer_props=DataTransferProps())
Definition H5Slice_traits_misc.hpp:256
Selection select(const HyperSlab &hyperslab) const
Select an hyperslab in the current Slice/Dataset.
Definition H5Slice_traits_misc.hpp:82
void read_raw(T *array, const DataType &dtype, const DataTransferProps &xfer_props=DataTransferProps()) const
Definition H5Slice_traits_misc.hpp:227
void write_raw(const T *buffer, const DataType &mem_datatype, const DataTransferProps &xfer_props=DataTransferProps())
Definition H5Slice_traits_misc.hpp:281
Definition H5_definitions.hpp:22
Definition H5Slice_traits.hpp:73