Theoretica
Scientific Computing
Loading...
Searching...
No Matches
hdf5.h
Go to the documentation of this file.
1
8
9#ifndef THEORETICA_IO_HDF5_H
10#define THEORETICA_IO_HDF5_H
11
12#include <hdf5.h>
13
14#include <string>
15#include <vector>
16#include <unordered_map>
17#include <stdexcept>
18#include <cstring>
19#include <sstream>
20#include <algorithm>
21
22#include "./error.h"
23#include "../core/constants.h"
24#include "../algebra/vec.h"
25#include "../algebra/mat.h"
26
27
28namespace theoretica {
29namespace io {
30
31
32 // Data structures and representations
33
35 enum class HDF5NodeType {
36
38 UNKNOWN,
39
41 GROUP,
42
45 };
46
47
50 struct hdf5_node {
51
53 std::string name;
54
57 std::string path;
58
61
63 std::vector<size_t> dimensions;
64
66 std::vector<std::string> attributes;
67
69 std::unordered_map<std::string, hdf5_node> children;
70
71
73 inline bool is_group() const noexcept {
74 return type == HDF5NodeType::GROUP;
75 }
76
77
79 inline bool is_dataset() const noexcept {
81 }
82
83
88 inline hdf5_node& operator[](const std::string& child_name) {
89 return children[child_name];
90 }
91
92
97 inline const hdf5_node& operator[](const std::string& child_name) const {
98 return children.at(child_name);
99 }
100 };
101
102
105 struct hdf5_handle {
106
110
111
116
117
120
121 if (id == H5I_INVALID_HID || !H5Iis_valid(id)) return;
122 switch (H5Iget_type(id)) {
123 case H5I_FILE: H5Fclose(id); break;
124 case H5I_GROUP: H5Gclose(id); break;
125 case H5I_DATASET: H5Dclose(id); break;
126 case H5I_DATATYPE: H5Tclose(id); break;
127 case H5I_DATASPACE: H5Sclose(id); break;
128 case H5I_ATTR: H5Aclose(id); break;
129 case H5I_GENPROP_LST: H5Pclose(id); break;
130 default: H5Idec_ref(id); break;
131 }
132 }
133
134
136 hdf5_handle(const hdf5_handle&) = delete;
137 hdf5_handle& operator=(const hdf5_handle&) = delete;
138
139
144
145
148
149 if (this != &other) {
150
151 if (id >= 0 && H5Iis_valid(id))
152 H5Idec_ref(id);
153
154 id = other.id;
156 }
157
158 return *this;
159 }
160
161
163 operator hid_t() const {
164 return id;
165 }
166 };
167
168
169 namespace _internal {
170
174 template<typename Type>
175 inline hid_t hdf5_type() {
176 TH_IO_ERROR("io::_internal::hdf5_type", "hdf5_type<Type>", IoError::FormatError);
177 return H5I_INVALID_HID;
178 }
179
180 template<> inline hid_t hdf5_type<double>() { return H5T_NATIVE_DOUBLE; }
181 template<> inline hid_t hdf5_type<float>() { return H5T_NATIVE_FLOAT; }
182 template<> inline hid_t hdf5_type<int>() { return H5T_NATIVE_INT; }
183 template<> inline hid_t hdf5_type<long>() { return H5T_NATIVE_LONG; }
184 template<> inline hid_t hdf5_type<unsigned int>() { return H5T_NATIVE_UINT; }
185
186
188 inline void suppress_errors() {
189 H5Eset_auto(H5E_DEFAULT, nullptr, nullptr);
190 }
191
192
196 inline void remove_link(const hdf5_handle& id, const std::string& path) {
197
199 if (H5Lexists(id, path.c_str(), H5P_DEFAULT) > 0) {
200 H5Ldelete(id, path.c_str(), H5P_DEFAULT);
201 }
202 }
203
204
205 // Attributes
206
207 // Count attributes callback for H5Aiterate2
208 static herr_t attribute_callback(hid_t id, const char* attr_name, const H5A_info_t* info_ptr, void* op_data) {
209
210 auto* attrs = (std::vector<std::string>*) op_data;
211 attrs->emplace_back(attr_name);
212 return 0;
213 }
214
215
216 // Load all attributes for a given object into a vector of strings
217 inline void load_attributes(hid_t obj_id, std::vector<std::string>& attributes) {
218 H5Aiterate2(obj_id, H5_INDEX_CRT_ORDER, H5_ITER_INC, NULL, attribute_callback, &attributes);
219 }
220
221
222 // Read and write attributes of various types, with error handling
223 template<typename Type>
224 inline void read_attribute(hid_t attr_id, Type& value) {
225
226 if (H5Aread(attr_id, hdf5_type<Type>(), &value) < 0) {
227 TH_IO_ERROR("io::read_attribute", "attribute_id", IoError::FormatError);
228 value = Type();
229 }
230 }
231
232
233 // Read a string attribute
234 inline void read_attribute(hid_t attr_id, std::string& value) {
235
236 hdf5_handle type_id = H5Aget_type(attr_id);
238 TH_IO_ERROR("io::read_attribute", "attribute_id", IoError::FormatError);
239 value = std::string();
240 return;
241 }
242
244 if (is_varstring > 0) {
245
246 char* buf = nullptr;
247 hdf5_handle mem_type = H5Tcopy(H5T_C_S1);
249
250 if (H5Aread(attr_id, mem_type, &buf) < 0) {
251 TH_IO_ERROR("io::read_attribute", "attribute_id", IoError::ReadError);
252 value = std::string();
253 return;
254 }
255
256 value = buf ? std::string(buf) : std::string();
257 if (buf)
259
260 } else if (is_varstring == 0) {
261
262 const size_t size = H5Tget_size(type_id);
263 std::vector<char> buf (size + 1, '\0');
264
265 hdf5_handle mem_type = H5Tcopy(H5T_C_S1);
266 H5Tset_size(mem_type, size);
267 if (H5Aread(attr_id, mem_type, buf.data()) < 0) {
268 TH_IO_ERROR("io::read_attribute", "attribute_id", IoError::ReadError);
269 value = std::string();
270 return;
271 }
272
273 value = std::string(buf.data());
274
275 } else {
276 TH_IO_ERROR("io::read_attribute", "attribute_id", IoError::FormatError);
277 value = std::string();
278 return;
279 }
280 }
281
282
283 // Write an attribute of a generic type, creating it if it doesn't exist
284 template<typename Type>
285 inline void write_attribute(hid_t obj_id, const std::string& name, const Type& value) {
286
287 hdf5_handle space_id = H5Screate(H5S_SCALAR);
288 hdf5_handle attr_id = H5Acreate2(obj_id, name.c_str(), hdf5_type<Type>(), space_id, H5P_DEFAULT, H5P_DEFAULT);
289
290 if (attr_id < 0) {
291 TH_IO_ERROR("io::write_attribute", "attribute_id", IoError::WriteError);
292 return;
293 }
294
295 if (H5Awrite(attr_id, hdf5_type<Type>(), &value) < 0)
296 TH_IO_ERROR("io::write_attribute", "attribute_id", IoError::WriteError);
297 }
298
299
300 // Write a string attribute
301 inline void write_attribute(hid_t obj_id, const std::string& name, const std::string& value) {
302
303 hdf5_handle space_id = H5Screate(H5S_SCALAR);
304 hdf5_handle type_id = H5Tcopy(H5T_C_S1);
305 H5Tset_size(type_id, value.empty() ? 1 : value.length());
307
308 hdf5_handle attr_id = H5Acreate2(obj_id, name.c_str(), type_id, space_id, H5P_DEFAULT, H5P_DEFAULT);
309 if (attr_id < 0) {
310 TH_IO_ERROR("io::write_attribute", "attribute_id", IoError::WriteError);
311 return;
312 }
313
314 if (H5Awrite(attr_id, type_id, value.c_str()) < 0) {
315 TH_IO_ERROR("io::write_attribute", "attribute_id", IoError::WriteError);
316 }
317 }
318
319
320 // Write a C-style string attribute
321 inline void write_attribute(hid_t obj_id, const std::string& name, const char* str) {
322 write_attribute(obj_id, name, std::string(str));
323 }
324
325
326 // Recursively iterates over group members to build the hdf5_node tree structure
327 inline herr_t iter_callback(hid_t id, const char* name, const H5L_info_t* info, void* opdata) {
328
329 hdf5_node* parent = (hdf5_node*) opdata;
330 hdf5_node child;
331 child.name = name;
332
333 // Build full child path
334 if (parent->path == "/")
335 child.path = "/" + std::string(name);
336 else
337 child.path = parent->path + "/" + std::string(name);
338
339 // Open object by absolute child path from current location
340 hdf5_handle obj_id = H5Oopen(id, name, H5P_DEFAULT);
341 if (obj_id < 0) {
343 parent->children.emplace(child.name, std::move(child));
344 return 0;
345 }
346
348 if (obj_type == H5I_GROUP) {
349
351
352 // Read group attributes
353 load_attributes(obj_id, child.attributes);
354
355 // Recursion by opening the group and iterating inside it
356 hdf5_handle gid = H5Gopen2(id, name, H5P_DEFAULT);
357 if (gid >= 0)
358 H5Literate(gid, H5_INDEX_NAME, H5_ITER_INC, NULL, iter_callback, &child);
359
360 } else if (obj_type == H5I_DATASET) {
361
363
364 // Dataset shape
365 hdf5_handle space_id = H5Dget_space(obj_id);
367 if (ndims > 0) {
368 std::vector<hsize_t> dims (static_cast<size_t>(ndims));
370 child.dimensions = dims;
371 }
372
373 // Dataset attributes
374 load_attributes(obj_id, child.attributes);
375
376 } else {
377
379 load_attributes(obj_id, child.attributes);
380 }
381
382 parent->children.emplace(child.name, std::move(child));
383 return 0;
384 }
385
386
387 // Helper function to recursively build a formatted
388 // tree string representation of the HDF5 structure
389 inline void build_tree(
390 std::ostringstream& oss, const hdf5_node& node,
391 const std::string& prefix, bool is_last, bool is_root) {
392
393 // Current node line
394 if (is_root)
395 oss << (node.name.empty() ? "/" : node.name);
396 else
397 oss << prefix << (is_last ? "└── " : "├── ") << node.name;
398
399
400 // Print dataset
401 if (node.is_dataset()) {
402
403 oss << " [";
404 for (size_t i = 0; i < node.dimensions.size(); ++i) {
405 oss << node.dimensions[i];
406 if (i + 1 < node.dimensions.size())
407 oss << ", ";
408 }
409 oss << "]";
410 }
411
412
413 // Show attribute count
414 if (!node.attributes.empty()) {
415 oss << " (@" << node.attributes.size() << " attrs)";
416 }
417 oss << "\n";
418
419
420 if (node.children.empty())
421 return;
422
423 // Print group children recursively
424 size_t i = 0;
425 const std::string child_prefix = is_root ? "" : (prefix + (is_last ? " " : "│ "));
426 for (auto& pair : node.children) {
427
428 const bool is_last_child = (node.children.size() == i + 1);
429 build_tree(
430 oss,
431 pair.second,
434 false
435 );
436
437 i++;
438 }
439 }
440 }
441
442
443 // Output formatting
444
449 inline std::string to_string(const hdf5_node& node) {
450 std::ostringstream oss;
451 _internal::build_tree(oss, node, "", true, true);
452 return oss.str();
453 }
454
455
461 inline std::ostream& operator<<(std::ostream& os, const hdf5_node& node) {
462 return os << to_string(node);
463 }
464
465
466 // Public API using handles, for continued operations on an open file or group
467
473 inline hdf5_handle hdf5_open(const std::string& filename, bool write = false) {
474
475 _internal::suppress_errors();
477
478 // Open the file with write permissions,
479 // creating it if it doesn't exist.
480 if (write) {
481
482 file_id = H5Fopen(filename.c_str(), H5F_ACC_RDWR, H5P_DEFAULT);
483
484 if (file_id < 0) {
485
487
488 if (file_id < 0) {
489 TH_IO_ERROR("io::hdf5_open", filename, IoError::WriteError);
490 return H5I_INVALID_HID;
491 }
492 }
493
494 } else {
495
496 // Open with read-only access.
497 file_id = H5Fopen(filename.c_str(), H5F_ACC_RDONLY, H5P_DEFAULT);
498
499 if (file_id < 0) {
500 TH_IO_ERROR("io::hdf5_open", filename, IoError::FileNotFound);
501 return H5I_INVALID_HID;
502 }
503 }
504
505 return hdf5_handle(file_id);
506 }
507
508
513 inline bool hdf5_is_valid(const hdf5_handle& handle) {
514 return H5Iis_valid(handle.id);
515 }
516
517
522 inline hdf5_node hdf5_load(const hdf5_handle& id) {
523
524 _internal::suppress_errors();
525
527 root.name = "/";
528 root.path = "/";
530 root.dimensions.clear();
531 root.attributes.clear();
532 root.children.clear();
533
534 _internal::load_attributes(id, root.attributes);
535 H5Literate(id, H5_INDEX_NAME, H5_ITER_INC, NULL, _internal::iter_callback, &root);
536
537 return root;
538 }
539
540
545 inline void hdf5_create_group(const hdf5_handle& id, const std::string& path) {
546
547 _internal::suppress_errors();
548
549 // Already exists, do nothing
550 if (H5Lexists(id, path.c_str(), H5P_DEFAULT) > 0)
551 return;
552
553 hid_t gid = H5Gcreate2(id, path.c_str(), H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
554 if (gid < 0) {
555 TH_IO_ERROR("io::hdf5_create_group", path, IoError::WriteError);
556 return;
557 }
558
560 }
561
562
567 inline void hdf5_delete_group(const hdf5_handle& id, const std::string& path) {
568
569 _internal::suppress_errors();
570
571 // Do nothing if it doesn't exist
572 if (H5Lexists(id, path.c_str(), H5P_DEFAULT) <= 0)
573 return;
574
575 if (H5Ldelete(id, path.c_str(), H5P_DEFAULT) < 0) {
576 TH_IO_ERROR("io::hdf5_delete_group", path, IoError::WriteError);
577 }
578 }
579
580
581 // Metadata handling
582
590 template<typename Type>
591 inline Type hdf5_read_attribute(const hdf5_handle& id, const std::string& path, const std::string& attr_name) {
592
593 _internal::suppress_errors();
594
595 hid_t obj_id = H5Oopen(id, path.c_str(), H5P_DEFAULT);
596 if (obj_id < 0) {
597 TH_IO_ERROR("io::read_attribute", path, IoError::FileNotFound);
598 return Type();
599 }
600
603 if (attr_id < 0) {
604 TH_IO_ERROR("io::read_attribute", path + "@" + attr_name, IoError::FileNotFound);
605 return Type();
606 }
608
609 Type value;
610 _internal::read_attribute(attr_id, value);
611 return value;
612 }
613
614
620 inline void hdf5_delete_attribute(const hdf5_handle& id, const std::string& path, const std::string& attr_name) {
621
622 _internal::suppress_errors();
623
624 hid_t obj_id = H5Oopen(id, path.c_str(), H5P_DEFAULT);
625 if (obj_id < 0) {
626 TH_IO_ERROR("io::delete_attribute", path, IoError::FileNotFound);
627 return;
628 }
629
631
632 if (H5Aexists(obj_id, attr_name.c_str()) > 0) {
633 if (H5Adelete(obj_id, attr_name.c_str()) < 0) {
634 TH_IO_ERROR("io::delete_attribute", path + "@" + attr_name, IoError::WriteError);
635 }
636 }
637 }
638
639
647 template<typename Type>
648 inline void hdf5_write_attribute(const hdf5_handle& id, const std::string& path, const std::string& attr_name, const Type& value) {
649
650 _internal::suppress_errors();
651 hid_t obj_id = H5Oopen(id, path.c_str(), H5P_DEFAULT);
652 if (obj_id < 0) {
653 TH_IO_ERROR("io::write_attribute", path, IoError::FileNotFound);
654 return;
655 }
656
658 if (H5Aexists(obj_id, attr_name.c_str()) > 0)
659 H5Adelete(obj_id, attr_name.c_str());
660
661 _internal::write_attribute(obj_id, attr_name, value);
662 }
663
664
665 // Dataset IO for vectors and matrices
666
667
674 template<typename Vector = vec<real>>
675 inline Vector& hdf5_read_vec(const hdf5_handle& id, const std::string& path, Vector& v) {
676
677 using Type = vector_element_t<Vector>;
678
679 _internal::suppress_errors();
680 hid_t data_id = H5Dopen2(id, path.c_str(), H5P_DEFAULT);
681 if (data_id < 0) {
682 TH_IO_ERROR("io::hdf5_read_vec", path, IoError::FileNotFound);
683 return algebra::vec_error(v);
684 }
686
688
689 // Dataset must be 1D to read into a vector
691 TH_IO_ERROR("io::hdf5_read_vec", path, IoError::FormatError);
692 return algebra::vec_error(v);
693 }
694
695 hsize_t dims[1];
697
698 v.resize(dims[0], nan());
699 if (v.size() != dims[0]) {
700 TH_IO_ERROR("io::hdf5_read_vec", path, IoError::FormatError);
701 return algebra::vec_error(v);
702 }
703
704 if (H5Dread(data_id, _internal::hdf5_type<Type>(), H5S_ALL, H5S_ALL, H5P_DEFAULT, v.data()) < 0) {
705 TH_IO_ERROR("io::hdf5_read_vec", path, IoError::ReadError);
706 return algebra::vec_error(v);
707 }
708
709 return v;
710 }
711
712
719 template<typename Vector = vec<real>>
720 inline Vector hdf5_read_vec(const hdf5_handle& id, const std::string& path) {
721 Vector v;
722 return hdf5_read_vec(id, path, v);
723 }
724
725
732 template<typename Vector>
733 inline void hdf5_write_vec(const hdf5_handle& id, const std::string& path, const Vector& v) {
734
735 using Type = vector_element_t<Vector>;
736 _internal::remove_link(id, path);
737
738 hsize_t dims[1] = { (hsize_t) v.size() };
740
741 hid_t data_id = H5Dcreate2(id, path.c_str(), _internal::hdf5_type<Type>(), h_space, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
742 if (data_id < 0) {
743 TH_IO_ERROR("io::hdf5_write_vec", path, IoError::WriteError);
744 return;
745 }
747
748 if (H5Dwrite(data_id, _internal::hdf5_type<Type>(), H5S_ALL, H5S_ALL, H5P_DEFAULT, v.data()) < 0) {
749 TH_IO_ERROR("io::hdf5_write_vec", path, IoError::WriteError);
750 }
751 }
752
753
758 inline void hdf5_delete_dataset(const hdf5_handle& id, const std::string& path) {
759
760 _internal::suppress_errors();
761
762 if (H5Lexists(id, path.c_str(), H5P_DEFAULT) > 0) {
763 if (H5Ldelete(id, path.c_str(), H5P_DEFAULT) < 0) {
764 TH_IO_ERROR("io::hdf5_delete_dataset", path, IoError::WriteError);
765 }
766 }
767 }
768
769
776 template<typename Matrix = mat<real>>
777 inline Matrix hdf5_read_mat(const hdf5_handle& id, const std::string& path, Matrix& m) {
778
779 using Type = matrix_element_t<Matrix>;
780 _internal::suppress_errors();
781
782 hid_t data_id = H5Dopen2(id, path.c_str(), H5P_DEFAULT);
783 if (data_id < 0) {
784 TH_IO_ERROR("io::hdf5_read_mat", path, IoError::FileNotFound);
785 return algebra::mat_error(m);
786 }
788
791 TH_IO_ERROR("io::hdf5_read_mat", path, IoError::FormatError);
792 return algebra::mat_error(m);
793 }
794
795 hsize_t dims[2];
797
798 m.resize(dims[0], dims[1]);
799 if (m.rows() != dims[0] || m.cols() != dims[1]) {
800 TH_IO_ERROR("io::hdf5_read_mat", path, IoError::FormatError);
801 return algebra::mat_error(m);
802 }
803
804 if (H5Dread(data_id, _internal::hdf5_type<Type>(), H5S_ALL, H5S_ALL, H5P_DEFAULT, m.data()) < 0) {
805 TH_IO_ERROR("io::hdf5_read_mat", path, IoError::ReadError);
806 return algebra::mat_error(m);
807 }
808
809 return m;
810 }
811
812
819 template<typename Matrix = mat<real>>
820 inline Matrix hdf5_read_mat(const hdf5_handle& id, const std::string& path) {
821 Matrix m;
822 return hdf5_read_mat(id, path, m);
823 }
824
825
832 template<typename Matrix>
833 inline void hdf5_write_mat(const hdf5_handle& id, const std::string& path, const Matrix& m) {
834
835 using Type = matrix_element_t<Matrix>;
836 _internal::remove_link(id, path);
837
838 hsize_t dims[2] = { hsize_t(m.rows()), hsize_t(m.cols()) };
840
841 hid_t data_id = H5Dcreate2(id, path.c_str(), _internal::hdf5_type<Type>(), h_space, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
842 if (data_id < 0) {
843 TH_IO_ERROR("io::hdf5_write_mat", "dataset_path", IoError::FormatError);
844 return;
845 }
846
848 if (H5Dwrite(data_id, _internal::hdf5_type<Type>(), H5S_ALL, H5S_ALL, H5P_DEFAULT, m.data()) < 0) {
849 TH_IO_ERROR("io::hdf5_write_mat", "dataset_path", IoError::WriteError);
850 }
851 }
852
853
856 class hdf5_file {
857 private:
858
860 std::string m_filename;
861
863 hdf5_handle m_file_id;
864
866 hdf5_node m_root;
867
868 public:
869
874 explicit hdf5_file(const std::string& filename, bool write = false) : m_filename(filename) {
875 m_file_id = hdf5_open(m_filename, write);
876
877 if (hdf5_is_valid(m_file_id))
878 refresh();
879 }
880
881
883 void refresh() {
884
885 // Always rebuild from scratch
886 m_root = hdf5_node{};
887 m_root.name = "/";
888 m_root.path = "/";
889 m_root.type = HDF5NodeType::GROUP;
890 m_root.dimensions.clear();
891 m_root.attributes.clear();
892 m_root.children.clear();
893
894 // Rebuild root from open file handle
895 m_root = hdf5_load(m_file_id);
896
897 // Keep root name as filename
898 m_root.name = m_filename;
899 }
900
901
903 void close() {
904
905 H5Fclose(m_file_id);
906 m_file_id = H5I_INVALID_HID;
907
908 // Clear cached structure
909 m_root = hdf5_node();
910 }
911
912
914 const std::string& filename() const {
915 return m_filename;
916 }
917
918
921 const hdf5_node& root() const {
922 return m_root;
923 }
924
925
928 hid_t id() const {
929 return m_file_id;
930 }
931
932
935 hdf5_node& operator[](const std::string& child_name) {
936 return m_root[child_name];
937 }
938
939
942 const hdf5_node& operator[](const std::string& child_name) const {
943 return m_root[child_name];
944 }
945
946
949 inline void create_group(const std::string& path) {
950 hdf5_create_group(m_file_id, path);
951 }
952
953
957 inline void create_group(const hdf5_node& parent, const std::string& child_name) {
958
959 const std::string path = (parent.path == "/")
960 ? ("/" + child_name)
961 : (parent.path + "/" + child_name);
962
963 hdf5_create_group(m_file_id, path);
964 }
965
966
969 inline void delete_group(const std::string& path) {
970 hdf5_delete_group(m_file_id, path);
971 }
972
973
976 inline void delete_group(const hdf5_node& node) {
977 hdf5_delete_group(m_file_id, node.path);
978 }
979
980
981 // Vector operations
982
989 template<typename Vector = vec<real>>
990 Vector read_vec(const std::string& path) const {
991 return hdf5_read_vec<Vector>(m_file_id, path);
992 }
993
994
1001 template<typename Vector = vec<real>>
1003 return hdf5_read_vec<Vector>(m_file_id, node.path);
1004 }
1005
1006
1012 template<typename Vector>
1013 void write_vec(const std::string& path, const Vector& v) {
1014 hdf5_write_vec(m_file_id, path, v);
1015 }
1016
1022 template<typename Vector>
1023 void write_vec(const hdf5_node& node, const Vector& v) {
1024 hdf5_write_vec(m_file_id, node.path, v);
1025 }
1026
1027
1028 // Matrix operations
1029
1035 template<typename Matrix = mat<real>>
1036 Matrix read_mat(const std::string& path) const {
1037 return hdf5_read_mat<Matrix>(m_file_id, path);
1038 }
1039
1040
1046 template<typename Matrix = mat<real>>
1048 return hdf5_read_mat<Matrix>(m_file_id, node.path);
1049 }
1050
1051
1057 template<typename Matrix>
1058 void write_mat(const std::string& path, const Matrix& m) {
1059 hdf5_write_mat(m_file_id, path, m);
1060 }
1061
1062
1068 template<typename Matrix>
1069 void write_mat(const hdf5_node& node, const Matrix& m) {
1070 hdf5_write_mat(m_file_id, node.path, m);
1071 }
1072
1073
1077 void delete_dataset(const std::string& path) {
1078 hdf5_delete_dataset(m_file_id, path);
1079 }
1080
1081
1086 hdf5_delete_dataset(m_file_id, node.path);
1087 }
1088
1089
1090 // Metadata operations
1091
1098 template<typename Type>
1099 Type read_attribute(const std::string& path, const std::string& attr_name) const {
1100 return hdf5_read_attribute<Type>(m_file_id, path, attr_name);
1101 }
1102
1103
1110 template<typename Type>
1111 Type read_attribute(const hdf5_node& node, const std::string& attr_name) const {
1112 return hdf5_read_attribute<Type>(m_file_id, node.path, attr_name);
1113 }
1114
1115
1122 template<typename Type>
1123 void write_attribute(const std::string& path, const std::string& attr_name, const Type& value) {
1124 hdf5_write_attribute<Type>(m_file_id, path, attr_name, value);
1125 }
1126
1127
1134 template<typename Type>
1135 void write_attribute(const hdf5_node& node, const std::string& attr_name, const Type& value) {
1136 hdf5_write_attribute<Type>(m_file_id, node.path, attr_name, value);
1137 }
1138
1139
1144 void delete_attribute(const std::string& path, const std::string& attr_name) {
1145 hdf5_delete_attribute(m_file_id, path, attr_name);
1146 }
1147
1148
1153 void delete_attribute(const hdf5_node& node, const std::string& attr_name) {
1154 hdf5_delete_attribute(m_file_id, node.path, attr_name);
1155 }
1156
1157
1162 inline std::string to_string() const {
1163 return io::to_string(m_root);
1164 }
1165
1166
1172 inline friend std::ostream& operator<<(std::ostream& os, const hdf5_file& file) {
1173 return os << file.to_string();
1174 }
1175 };
1176
1177}}
1178
1179#endif
High-level interface for managing an HDF5 file, its structure, and operations.
Definition hdf5.h:856
hid_t id() const
Get the raw ID of the HDF5 file handle for direct API access.
Definition hdf5.h:928
Matrix read_mat(const std::string &path) const
Read a 2D dataset array into a matrix from the given path.
Definition hdf5.h:1036
void write_vec(const std::string &path, const Vector &v)
Writes a 1D vector to an HDF5 dataset at the given path, overwriting if it exists.
Definition hdf5.h:1013
Type read_attribute(const hdf5_node &node, const std::string &attr_name) const
Read metadata (attribute) attached to a specific node by reference.
Definition hdf5.h:1111
void write_mat(const hdf5_node &node, const Matrix &m)
Write a 2D matrix to an HDF5 dataset at the given node, overwriting if it exists.
Definition hdf5.h:1069
void delete_attribute(const std::string &path, const std::string &attr_name)
Delete an attribute attached to a specific node by path, if it exists.
Definition hdf5.h:1144
void delete_group(const hdf5_node &node)
Delete a group referenced by a node.
Definition hdf5.h:976
const hdf5_node & root() const
Get the base node containing the entire loaded file structure.
Definition hdf5.h:921
void delete_dataset(const std::string &path)
Delete a dataset at the given path if it exists.
Definition hdf5.h:1077
void create_group(const std::string &path)
Create a new group inside the file.
Definition hdf5.h:949
void create_group(const hdf5_node &parent, const std::string &child_name)
Create a new group using a parent node and a child name.
Definition hdf5.h:957
void write_attribute(const std::string &path, const std::string &attr_name, const Type &value)
Write or overwrite metadata (attribute) attached to a specific node by path.
Definition hdf5.h:1123
void write_attribute(const hdf5_node &node, const std::string &attr_name, const Type &value)
Write or overwrite metadata (attribute) attached to a specific node by reference.
Definition hdf5.h:1135
void write_mat(const std::string &path, const Matrix &m)
Write a 2D matrix to an HDF5 dataset at the given path, overwriting if it exists.
Definition hdf5.h:1058
void delete_dataset(const hdf5_node &node)
Delete a dataset at the given node if it exists.
Definition hdf5.h:1085
void delete_group(const std::string &path)
Delete a group inside the file.
Definition hdf5.h:969
void close()
Close the HDF5 file immediately.
Definition hdf5.h:903
Vector read_vec(const hdf5_node &node) const
Loads a 1D dataset array into a vector from the given node inside the HDF5 file.
Definition hdf5.h:1002
hdf5_file(const std::string &filename, bool write=false)
Open the file persistently.
Definition hdf5.h:874
void write_vec(const hdf5_node &node, const Vector &v)
Writes a 1D vector to the HDF5 dataset represented by the given node, overwriting if it exists.
Definition hdf5.h:1023
Matrix read_mat(const hdf5_node &node) const
Read a 2D dataset array into a matrix from the given node.
Definition hdf5.h:1047
hdf5_node & operator[](const std::string &child_name)
Access children of the root node by reference using dictionary-like syntax.
Definition hdf5.h:935
void refresh()
Refresh the cached tree structure from disk.
Definition hdf5.h:883
Type read_attribute(const std::string &path, const std::string &attr_name) const
Read metadata (attribute) attached to a specific node by path.
Definition hdf5.h:1099
const std::string & filename() const
Get the name of the file.
Definition hdf5.h:914
Vector read_vec(const std::string &path) const
Loads a 1D dataset array into a vector from the given path inside the HDF5 file.
Definition hdf5.h:990
std::string to_string() const
Converts an entire HDF5 file structure into a formatted tree string.
Definition hdf5.h:1162
friend std::ostream & operator<<(std::ostream &os, const hdf5_file &file)
Print the HDF5 file structure to a stream in a formatted tree representation.
Definition hdf5.h:1172
const hdf5_node & operator[](const std::string &child_name) const
Access children of the root node by const reference using dictionary-like syntax.
Definition hdf5.h:942
void delete_attribute(const hdf5_node &node, const std::string &attr_name)
Delete an attribute attached to a specific node by path, if it exists.
Definition hdf5.h:1153
HDF5 file IO.
void remove_link(const hdf5_handle &id, const std::string &path)
Removes a link (dataset or group) if it exists at the given path.
Definition hdf5.h:196
void suppress_errors()
Suppress HDF5 error messages by setting a custom dummy error handler.
Definition hdf5.h:188
hid_t hdf5_type()
Helper template to map C++ types to HDF5 Native Types.
Definition hdf5.h:175
Error handling for IO operations.
Matrix & mat_error(Matrix &m)
Overwrite the given matrix with the error matrix with NaN values on the diagonal and zeroes everywher...
Definition algebra.h:40
Vector & vec_error(Vector &v)
Overwrite the given vector with the error vector with NaN values.
Definition algebra.h:58
hdf5_node hdf5_load(const hdf5_handle &id)
Recursively loads the structure of an active HDF5 group or file.
Definition hdf5.h:522
void hdf5_create_group(const hdf5_handle &id, const std::string &path)
Create a new group at the given path under an already open HDF5 location.
Definition hdf5.h:545
void hdf5_delete_dataset(const hdf5_handle &id, const std::string &path)
Deletes a dataset at the given path if it exists.
Definition hdf5.h:758
void hdf5_delete_group(const hdf5_handle &id, const std::string &path)
Delete a group (link) at the given path under an already open HDF5 location.
Definition hdf5.h:567
std::ostream & operator<<(std::ostream &os, const hdf5_node &node)
Prints the HDF5 file tree to a stream.
Definition hdf5.h:461
bool hdf5_is_valid(const hdf5_handle &handle)
Check whether a given HDF5 handle is valid.
Definition hdf5.h:513
Matrix hdf5_read_mat(const hdf5_handle &id, const std::string &path, Matrix &m)
Loads a 2D dataset array into a matrix.
Definition hdf5.h:777
@ FileNotFound
File or directory not found.
@ WriteError
Error occurred while writing to the file or stream.
@ FormatError
The file format is invalid or the data is corrupted.
@ ReadError
Error occurred while reading from the file or stream.
Type hdf5_read_attribute(const hdf5_handle &id, const std::string &path, const std::string &attr_name)
Reads an attribute attached to a specific node.
Definition hdf5.h:591
void hdf5_write_vec(const hdf5_handle &id, const std::string &path, const Vector &v)
Writes a 1D vector to an HDF5 dataset, overwriting if it exists.
Definition hdf5.h:733
void hdf5_delete_attribute(const hdf5_handle &id, const std::string &path, const std::string &attr_name)
Deletes an attribute attached to a specific node, if it exists.
Definition hdf5.h:620
Vector & hdf5_read_vec(const hdf5_handle &id, const std::string &path, Vector &v)
Loads a 1D dataset array into a vector.
Definition hdf5.h:675
HDF5NodeType
Type of node inside an HDF5 file.
Definition hdf5.h:35
@ DATASET
A node containing multi-dimensional array data.
@ UNKNOWN
An unsupported node type.
@ GROUP
A structural grouping containing other groups or datasets.
void hdf5_write_mat(const hdf5_handle &id, const std::string &path, const Matrix &m)
Writes a 2D matrix to an HDF5 dataset, overwriting if it exists.
Definition hdf5.h:833
void hdf5_write_attribute(const hdf5_handle &id, const std::string &path, const std::string &attr_name, const Type &value)
Writes or overwrites an attribute attached to a specific node.
Definition hdf5.h:648
hdf5_handle hdf5_open(const std::string &filename, bool write=false)
Open an HDF5 file with the given filename, returning a file handle.
Definition hdf5.h:473
Main namespace of the library which contains all functions and objects.
Definition algebra.h:27
Vector make_error()
Create a vector representing an error state, with all NaN values.
Definition algebra.h:103
TH_CONSTEXPR real nan()
Return a quiet NaN number in floating point representation.
Definition error.h:74
std::string to_string(MathError err)
Convert a MathError class enum to a string description.
Definition error.h:68
real root(real x, int n)
Compute the n-th root of x.
Definition real_analysis.h:812
RAII wrapper for managing HDF5 C-style handles, allowing direct API usage.
Definition hdf5.h:105
hid_t id
Underlying ID obtained by opening a file or group.
Definition hdf5.h:109
hdf5_handle(const hdf5_handle &)=delete
Prevent copying.
hdf5_handle(hdf5_handle &&other) noexcept
Move constructor.
Definition hdf5.h:141
hdf5_handle & operator=(hdf5_handle &&other) noexcept
Move assignment.
Definition hdf5.h:147
~hdf5_handle()
Safely decrements reference count and closes the handle if valid.
Definition hdf5.h:119
hdf5_handle(hid_t id=H5I_INVALID_HID)
Initialize the handle with a given HDF5 ID, for example obtained by opening a file or group,...
Definition hdf5.h:115
Represents a single node (group or dataset) in the HDF5 file hierarchy.
Definition hdf5.h:50
std::unordered_map< std::string, hdf5_node > children
For groups, the child nodes indexed by their local name.
Definition hdf5.h:69
std::vector< std::string > attributes
List of metadata attribute names attached to the node.
Definition hdf5.h:66
bool is_group() const noexcept
Check whether the node is a group.
Definition hdf5.h:73
std::vector< size_t > dimensions
For datasets, the dimensions of the array.
Definition hdf5.h:63
bool is_dataset() const noexcept
Check whether the node is a dataset.
Definition hdf5.h:79
HDF5NodeType type
The type of the node.
Definition hdf5.h:60
const hdf5_node & operator[](const std::string &child_name) const
Read-only dictionary-like access to child nodes by name.
Definition hdf5.h:97
std::string name
The local name of the node (e.g. "my_dataset")
Definition hdf5.h:53
hdf5_node & operator[](const std::string &child_name)
Dictionary-like access to child nodes by name.
Definition hdf5.h:88
std::string path
The absolute internal path to the node (e.g.
Definition hdf5.h:57