HighFive 2.10.0
HighFive - Header-only C++ HDF5 interface
Loading...
Searching...
No Matches
H5Easy_scalar.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
14namespace H5Easy {
15
16namespace detail {
17
18/*
19Base template for partial specialization: the fallback if specialized templates don't match.
20Used e.g. for scalars.
21*/
22template <typename T, typename = void>
23struct io_impl {
24 inline static DataSet dump(File& file,
25 const std::string& path,
26 const T& data,
27 const DumpOptions& options) {
28 DataSet dataset = initScalarDataset(file, path, data, options);
29 dataset.write(data);
30 if (options.flush()) {
31 file.flush();
32 }
33 return dataset;
34 }
35
36 inline static T load(const File& file, const std::string& path) {
37 DataSet dataset = file.getDataSet(path);
38 T data;
39 dataset.read(data);
40 return data;
41 }
42
43 inline static Attribute dumpAttribute(File& file,
44 const std::string& path,
45 const std::string& key,
46 const T& data,
47 const DumpOptions& options) {
48 Attribute attribute = initScalarAttribute(file, path, key, data, options);
49 attribute.write(data);
50 if (options.flush()) {
51 file.flush();
52 }
53 return attribute;
54 }
55
56 inline static T loadAttribute(const File& file,
57 const std::string& path,
58 const std::string& key) {
59 DataSet dataset = file.getDataSet(path);
60 Attribute attribute = dataset.getAttribute(key);
61 T data;
62 attribute.read(data);
63 return data;
64 }
65
66 inline static DataSet dump_extend(File& file,
67 const std::string& path,
68 const T& data,
69 const std::vector<size_t>& idx,
70 const DumpOptions& options) {
71 std::vector<size_t> ones(idx.size(), 1);
72
73 if (file.exist(path)) {
74 DataSet dataset = file.getDataSet(path);
75 std::vector<size_t> dims = dataset.getDimensions();
76 std::vector<size_t> shape = dims;
77 if (dims.size() != idx.size()) {
78 throw detail::error(
79 file,
80 path,
81 "H5Easy::dump: Dimension of the index and the existing field do not match");
82 }
83 for (size_t i = 0; i < dims.size(); ++i) {
84 shape[i] = std::max(dims[i], idx[i] + 1);
85 }
86 if (shape != dims) {
87 dataset.resize(shape);
88 }
89 dataset.select(idx, ones).write(data);
90 if (options.flush()) {
91 file.flush();
92 }
93 return dataset;
94 }
95
96 std::vector<size_t> shape = idx;
97 const size_t unlim = DataSpace::UNLIMITED;
98 std::vector<size_t> unlim_shape(idx.size(), unlim);
99 std::vector<hsize_t> chunks(idx.size(), 10);
100 if (options.isChunked()) {
101 chunks = options.getChunkSize();
102 if (chunks.size() != idx.size()) {
103 throw error(file, path, "H5Easy::dump: Incorrect dimension ChunkSize");
104 }
105 }
106 for (size_t& i: shape) {
107 i++;
108 }
109 DataSpace dataspace = DataSpace(shape, unlim_shape);
110 DataSetCreateProps props;
111 props.add(Chunking(chunks));
112 DataSet dataset = file.createDataSet(path, dataspace, AtomicType<T>(), props, {}, true);
113 dataset.select(idx, ones).write(data);
114 if (options.flush()) {
115 file.flush();
116 }
117 return dataset;
118 }
119
120 inline static T load_part(const File& file,
121 const std::string& path,
122 const std::vector<size_t>& idx) {
123 std::vector<size_t> ones(idx.size(), 1);
124 DataSet dataset = file.getDataSet(path);
125 T data;
126 dataset.select(idx, ones).read(data);
127 return data;
128 }
129};
130
131} // namespace detail
132} // namespace H5Easy
static const size_t UNLIMITED
Magic value to specify that a DataSpace can grow without limit.
Definition H5DataSpace.hpp:49
PropertyList< PropertyType::DATASET_CREATE > DataSetCreateProps
Definition H5PropertyList.hpp:201
Read/dump DataSets or Attribute using a minimalistic syntax. To this end, the functions are templated...
Definition H5Easy.hpp:59
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