GNU Octave 10.1.0
A high-level interpreted language, primarily intended for numerical computations, mostly compatible with Matlab
 
Loading...
Searching...
No Matches
ls-hdf5.cc
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////
2//
3// Copyright (C) 1996-2025 The Octave Project Developers
4//
5// See the file COPYRIGHT.md in the top-level directory of this
6// distribution or <https://octave.org/copyright/>.
7//
8// This file is part of Octave.
9//
10// Octave is free software: you can redistribute it and/or modify it
11// under the terms of the GNU General Public License as published by
12// the Free Software Foundation, either version 3 of the License, or
13// (at your option) any later version.
14//
15// Octave is distributed in the hope that it will be useful, but
16// WITHOUT ANY WARRANTY; without even the implied warranty of
17// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18// GNU General Public License for more details.
19//
20// You should have received a copy of the GNU General Public License
21// along with Octave; see the file COPYING. If not, see
22// <https://www.gnu.org/licenses/>.
23//
24////////////////////////////////////////////////////////////////////////
25
26#if defined (HAVE_CONFIG_H)
27# include "config.h"
28#endif
29
30#if defined (HAVE_HDF5)
31
32#include <cctype>
33
34#include <iomanip>
35#include <istream>
36#include <limits>
37#include <ostream>
38#include <string>
39#include <vector>
40
41#include "byte-swap.h"
42#include "data-conv.h"
43#include "file-ops.h"
44#include "glob-match.h"
45#include "lo-mappers.h"
46#include "mach-info.h"
47#include "oct-env.h"
48#include "oct-time.h"
49#include "quit.h"
50#include "str-vec.h"
51#include "oct-locbuf.h"
52
53#include "Cell.h"
54#include "defun.h"
55#include "error.h"
56#include "errwarn.h"
57#include "interpreter.h"
58#include "interpreter-private.h"
59#include "load-save.h"
60#include "oct-hdf5.h"
61#include "ovl.h"
62#include "oct-map.h"
63#include "ov-cell.h"
64#include "pager.h"
65#include "sysdep.h"
66#include "unwind-prot.h"
67#include "utils.h"
68#include "variables.h"
69#include "version.h"
70#include "dMatrix.h"
71#include "ov-lazy-idx.h"
72
73#include "ls-utils.h"
74#include "ls-hdf5.h"
75
76#if defined (HAVE_HDF5)
77
78static hid_t
79check_hdf5_id_value (octave_hdf5_id id, const char *who)
80{
81 if (id > std::numeric_limits<hid_t>::max ())
82 error ("%s: internal error: ID too large for hid_t", who);
83
84 return static_cast<hid_t> (id);
85}
86
87#else
88
89OCTAVE_NORETURN static void
90error_unexpected (const char *name)
91{
92 error ("unexpected call to %s when HAVE_HDF5 is not defined - please report this bug", name);
93}
94
95#endif
96
97hdf5_fstreambase::hdf5_fstreambase (const char *name, int mode, int /* prot */)
98 : file_id (-1), current_item (-1)
99{
100#if defined (HAVE_HDF5)
101 open_create (name, mode);
102
103 current_item = 0;
104
105#else
106 err_disabled_feature ("hdf5_fstreambase", "HDF5");
107#endif
108}
109
110void
112{
113#if defined (HAVE_HDF5)
114
115 if (file_id >= 0)
116 {
117 if (H5Fclose (file_id) < 0)
118 std::ios::setstate (std::ios::badbit);
119 file_id = -1;
120 }
121
122#else
123
124 error_unexpected ("hdf5_fstreambase::close");
125
126#endif
127}
128
129void
130hdf5_fstreambase::open (const char *name, int mode, int)
131{
132#if defined (HAVE_HDF5)
133
134 clear ();
135
136 open_create (name, mode);
137
138 current_item = 0;
139
140#else
141
142 error_unexpected ("hdf5_fstreambase::open");
143
144#endif
145}
146
147void
148hdf5_fstreambase::open_create (const char *name, int mode)
149{
150#if defined (HAVE_HDF5)
151 // Open the HDF5 file NAME. If it does not exist, create the file.
152
153# if defined (HAVE_HDF5_UTF8)
154 const char *fname = name;
155# else
156 std::string fname_str (name);
157 std::string ascii_fname_str = octave::sys::get_ASCII_filename (fname_str);
158 const char *fname = ascii_fname_str.c_str ();
159# endif
160
161 if (mode & std::ios::in)
162 file_id = H5Fopen (fname, H5F_ACC_RDONLY, octave_H5P_DEFAULT);
163 else if (mode & std::ios::out)
164 {
165 if (mode & std::ios::app && H5Fis_hdf5 (fname) > 0)
166 file_id = H5Fopen (fname, H5F_ACC_RDWR, octave_H5P_DEFAULT);
167 else
168# if defined (HAVE_HDF5_UTF8)
169 file_id = H5Fcreate (fname, H5F_ACC_TRUNC, octave_H5P_DEFAULT,
171# else
172 {
173 // Check whether file already exists
174 std::string abs_ascii_fname
175 = octave::sys::canonicalize_file_name (ascii_fname_str);
176 if (! abs_ascii_fname.empty ())
177 {
178 // Use the existing file
179 file_id = H5Fcreate (fname, H5F_ACC_TRUNC,
181 if (file_id < 0)
182 std::ios::setstate (std::ios::badbit);
183
184 return;
185 }
186
187 // Check whether filename contains non-ASCII (UTF-8) characters.
188 std::string::const_iterator first_non_ASCII
189 = std::find_if (fname_str.begin (), fname_str.end (),
190 [](char c) { return (c < 0 || c >= 128); });
191 if (first_non_ASCII == fname_str.end ())
192 {
193 // No non-ASCII characters
194 file_id = H5Fcreate (name, H5F_ACC_TRUNC, octave_H5P_DEFAULT,
196 if (file_id < 0)
197 std::ios::setstate (std::ios::badbit);
198
199 return;
200 }
201
202 // Create file in temp folder
203 std::string tmp_name = octave::sys::tempnam ("", "oct-");
204 octave_hdf5_id hdf5_fid = H5Fcreate (tmp_name.c_str (), H5F_ACC_TRUNC,
207 if (hdf5_fid < 0)
208 {
209 file_id = -1;
210 std::ios::setstate (std::ios::badbit);
211 return;
212 }
213
214 // Close file
215 H5Fclose (hdf5_fid);
216
217 // Move temporary file to final destination
218 std::string msg;
219 int res = octave::sys::rename (tmp_name, name, msg);
220 if (res < 0)
221 {
222 std::ios::setstate (std::ios::badbit);
223 file_id = -1;
224 return;
225 }
226
227 // Open file at final location
228 ascii_fname_str = octave::sys::get_ASCII_filename (fname_str);
229 ascii_fname = ascii_fname_str.c_str ();
230 file_id = H5Fopen (ascii_fname, H5F_ACC_RDWR, octave_H5P_DEFAULT);
231 }
232# endif
233 }
234 if (file_id < 0)
235 std::ios::setstate (std::ios::badbit);
236
237 return;
238
239#else
240
241 error_unexpected ("hdf5_fstreambase::open_create");
242
243#endif
244}
245
246static std::string
247make_valid_identifier (const std::string& nm)
248{
249 std::string retval;
250
251 std::size_t nm_len = nm.length ();
252
253 if (nm_len > 0)
254 {
255 if (! isalpha (nm[0]))
256 retval += '_';
257
258 for (std::size_t i = 0; i < nm_len; i++)
259 {
260 char c = nm[i];
261 retval += (isalnum (c) || c == '_') ? c : '_';
262 }
263 }
264
265 return retval;
266}
267
268// Given two compound types t1 and t2, determine whether they
269// are compatible for reading/writing. This function only
270// works for non-nested types composed of simple elements (ints, floats...),
271// which is all we need it for
272
273bool
275{
276#if defined (HAVE_HDF5)
277
278 int n;
279 if ((n = H5Tget_nmembers (t1)) != H5Tget_nmembers (t2))
280 return false;
281
282 for (int i = 0; i < n; ++i)
283 {
284 hid_t mt1 = H5Tget_member_type (t1, i);
285 hid_t mt2 = H5Tget_member_type (t2, i);
286
287 if (H5Tget_class (mt1) != H5Tget_class (mt2))
288 return false;
289
290 H5Tclose (mt2);
291 H5Tclose (mt1);
292 }
293
294 return true;
295
296#else
297 err_disabled_feature ("hdf5_types_compatible", "HDF5");
298#endif
299}
300
301// Return true if loc_id has the attribute named attr_name, and false
302// otherwise.
303
304bool
305hdf5_check_attr (octave_hdf5_id loc_id, const char *attr_name)
306{
307#if defined (HAVE_HDF5)
308
309 bool retval = false;
310
311 // we have to pull some shenanigans here to make sure
312 // HDF5 doesn't print out all sorts of error messages if we
313 // call H5Aopen for a non-existing attribute
314
315 H5E_auto_t err_fcn;
316 void *err_fcn_data;
317
318 // turn off error reporting temporarily, but save the error
319 // reporting function:
320
321#if defined (HAVE_HDF5_18)
322 H5Eget_auto (octave_H5E_DEFAULT, &err_fcn, &err_fcn_data);
323 H5Eset_auto (octave_H5E_DEFAULT, nullptr, nullptr);
324#else
325 H5Eget_auto (&err_fcn, &err_fcn_data);
326 H5Eset_auto (nullptr, nullptr);
327#endif
328
329 hid_t attr_id = H5Aopen_name (loc_id, attr_name);
330
331 if (attr_id >= 0)
332 {
333 // successful
334 retval = true;
335 H5Aclose (attr_id);
336 }
337
338 // restore error reporting:
339#if defined (HAVE_HDF5_18)
340 H5Eset_auto (octave_H5E_DEFAULT, err_fcn, err_fcn_data);
341#else
342 H5Eset_auto (err_fcn, err_fcn_data);
343#endif
344 return retval;
345
346#else
347 err_disabled_feature ("hdf5_check_attr", "HDF5");
348#endif
349}
350
351bool
353 const char *attr_name, void *buf)
354{
355#if defined (HAVE_HDF5)
356
357 bool retval = false;
358
359 // we have to pull some shenanigans here to make sure
360 // HDF5 doesn't print out all sorts of error messages if we
361 // call H5Aopen for a non-existing attribute
362
363 H5E_auto_t err_fcn;
364 void *err_fcn_data;
365
366 // turn off error reporting temporarily, but save the error
367 // reporting function:
368
369#if defined (HAVE_HDF5_18)
370 H5Eget_auto (octave_H5E_DEFAULT, &err_fcn, &err_fcn_data);
371 H5Eset_auto (octave_H5E_DEFAULT, nullptr, nullptr);
372#else
373 H5Eget_auto (&err_fcn, &err_fcn_data);
374 H5Eset_auto (nullptr, nullptr);
375#endif
376
377 hid_t attr_id = H5Aopen_name (loc_id, attr_name);
378
379 if (attr_id >= 0)
380 {
381 hid_t space_id = H5Aget_space (attr_id);
382
383 hsize_t rank = H5Sget_simple_extent_ndims (space_id);
384
385 if (rank == 0)
386 retval = H5Aread (attr_id, type_id, buf) >= 0;
387 H5Aclose (attr_id);
388 }
389
390 // restore error reporting:
391#if defined (HAVE_HDF5_18)
392 H5Eset_auto (octave_H5E_DEFAULT, err_fcn, err_fcn_data);
393#else
394 H5Eset_auto (err_fcn, err_fcn_data);
395#endif
396 return retval;
397
398#else
399 err_disabled_feature ("hdf5_get_scalar_attr", "HDF5");
400#endif
401}
402
403// The following subroutines creates an HDF5 representations of the way
404// we will store Octave complex types (pairs of floating-point numbers).
405// NUM_TYPE is the HDF5 numeric type to use for storage (e.g.
406// H5T_NATIVE_DOUBLE to save as 'double'). Note that any necessary
407// conversions are handled automatically by HDF5.
408
411{
412#if defined (HAVE_HDF5)
413
414 hid_t type_id = H5Tcreate (H5T_COMPOUND, sizeof (double) * 2);
415
416 H5Tinsert (type_id, "real", 0 * sizeof (double), num_type);
417 H5Tinsert (type_id, "imag", 1 * sizeof (double), num_type);
418
419 return type_id;
420
421#else
422 err_disabled_feature ("hdf5_make_complex_type", "HDF5");
423#endif
424}
425
426#if defined (HAVE_HDF5)
427
428// The following subroutine creates an HDF5 representation of the way
429// we will store Octave range types (triplets of floating-point numbers).
430// NUM_TYPE is the HDF5 numeric type to use for storage
431// (e.g., H5T_NATIVE_DOUBLE to save as 'double').
432// Note that any necessary conversions are handled automatically by HDF5.
433
434static hid_t
435hdf5_make_range_type (hid_t num_type)
436{
437 hid_t type_id = H5Tcreate (H5T_COMPOUND, sizeof (double) * 3);
438
439 H5Tinsert (type_id, "base", 0 * sizeof (double), num_type);
440 H5Tinsert (type_id, "limit", 1 * sizeof (double), num_type);
441 H5Tinsert (type_id, "increment", 2 * sizeof (double), num_type);
442
443 return type_id;
444}
445
446static herr_t
447load_inline_fcn (hid_t loc_id, const char *name, octave_value& retval)
448{
449#if defined (HAVE_HDF5)
450
451 hid_t group_hid, data_hid, space_hid, type_hid, type_class_hid, st_id;
452 hsize_t rank;
453 int slen;
454
455#if defined (HAVE_HDF5_18)
456 group_hid = H5Gopen (loc_id, name, octave_H5P_DEFAULT);
457#else
458 group_hid = H5Gopen (loc_id, name);
459#endif
460 if (group_hid < 0) return -1;
461
462#if defined (HAVE_HDF5_18)
463 data_hid = H5Dopen (group_hid, "args", octave_H5P_DEFAULT);
464#else
465 data_hid = H5Dopen (group_hid, "args");
466#endif
467 space_hid = H5Dget_space (data_hid);
468 rank = H5Sget_simple_extent_ndims (space_hid);
469
470 if (rank != 2)
471 {
472 H5Dclose (data_hid);
473 H5Sclose (space_hid);
474 H5Gclose (group_hid);
475 return -1;
476 }
477
478 OCTAVE_LOCAL_BUFFER (hsize_t, hdims, rank);
479 OCTAVE_LOCAL_BUFFER (hsize_t, maxdims, rank);
480
481 H5Sget_simple_extent_dims (space_hid, hdims, maxdims);
482
483 octave_value_list args (hdims[1]+1);
484
485 OCTAVE_LOCAL_BUFFER (char, s1, hdims[0] * hdims[1]);
486
487 if (H5Dread (data_hid, H5T_NATIVE_UCHAR, octave_H5S_ALL, octave_H5S_ALL,
488 octave_H5P_DEFAULT, s1) < 0)
489 {
490 H5Dclose (data_hid);
491 H5Sclose (space_hid);
492 H5Gclose (group_hid);
493 return -1;
494 }
495
496 H5Dclose (data_hid);
497 H5Sclose (space_hid);
498
499 for (std::size_t i = 0; i < hdims[1]; i++)
500 args(i+1) = std::string (s1 + i*hdims[0]);
501
502#if defined (HAVE_HDF5_18)
503 data_hid = H5Dopen (group_hid, "nm", octave_H5P_DEFAULT);
504#else
505 data_hid = H5Dopen (group_hid, "nm");
506#endif
507
508 if (data_hid < 0)
509 {
510 H5Gclose (group_hid);
511 return -1;
512 }
513
514 type_hid = H5Dget_type (data_hid);
515 type_class_hid = H5Tget_class (type_hid);
516
517 if (type_class_hid != H5T_STRING)
518 {
519 H5Tclose (type_hid);
520 H5Dclose (data_hid);
521 H5Gclose (group_hid);
522 return -1;
523 }
524
525 space_hid = H5Dget_space (data_hid);
526 rank = H5Sget_simple_extent_ndims (space_hid);
527
528 if (rank != 0)
529 {
530 H5Sclose (space_hid);
531 H5Tclose (type_hid);
532 H5Dclose (data_hid);
533 H5Gclose (group_hid);
534 return -1;
535 }
536
537 slen = H5Tget_size (type_hid);
538 if (slen < 0)
539 {
540 H5Sclose (space_hid);
541 H5Tclose (type_hid);
542 H5Dclose (data_hid);
543 H5Gclose (group_hid);
544 return -1;
545 }
546
547 OCTAVE_LOCAL_BUFFER (char, nm_tmp, slen);
548
549 // create datatype for (null-terminated) string to read into:
550 st_id = H5Tcopy (H5T_C_S1);
551 H5Tset_size (st_id, slen);
552
553 if (H5Dread (data_hid, st_id, octave_H5S_ALL, octave_H5S_ALL,
554 octave_H5P_DEFAULT, nm_tmp) < 0)
555 {
556 H5Sclose (space_hid);
557 H5Tclose (type_hid);
558 H5Gclose (group_hid);
559 return -1;
560 }
561 H5Tclose (st_id);
562 H5Dclose (data_hid);
563 // NAME is obsolete and unused.
564 // std::string name (nm_tmp);
565
566#if defined (HAVE_HDF5_18)
567 data_hid = H5Dopen (group_hid, "iftext", octave_H5P_DEFAULT);
568#else
569 data_hid = H5Dopen (group_hid, "iftext");
570#endif
571
572 if (data_hid < 0)
573 {
574 H5Gclose (group_hid);
575 return -1;
576 }
577
578 type_hid = H5Dget_type (data_hid);
579 type_class_hid = H5Tget_class (type_hid);
580
581 if (type_class_hid != H5T_STRING)
582 {
583 H5Tclose (type_hid);
584 H5Dclose (data_hid);
585 H5Gclose (group_hid);
586 return -1;
587 }
588
589 space_hid = H5Dget_space (data_hid);
590 rank = H5Sget_simple_extent_ndims (space_hid);
591
592 if (rank != 0)
593 {
594 H5Sclose (space_hid);
595 H5Tclose (type_hid);
596 H5Dclose (data_hid);
597 H5Gclose (group_hid);
598 return -1;
599 }
600
601 slen = H5Tget_size (type_hid);
602 if (slen < 0)
603 {
604 H5Sclose (space_hid);
605 H5Tclose (type_hid);
606 H5Dclose (data_hid);
607 H5Gclose (group_hid);
608 return -1;
609 }
610
611 OCTAVE_LOCAL_BUFFER (char, iftext_tmp, slen);
612
613 // create datatype for (null-terminated) string to read into:
614 st_id = H5Tcopy (H5T_C_S1);
615 H5Tset_size (st_id, slen);
616
617 if (H5Dread (data_hid, st_id, octave_H5S_ALL, octave_H5S_ALL,
618 octave_H5P_DEFAULT, iftext_tmp) < 0)
619 {
620 H5Sclose (space_hid);
621 H5Tclose (type_hid);
622 H5Gclose (group_hid);
623 return -1;
624 }
625 H5Tclose (st_id);
626 H5Dclose (data_hid);
627
628 args(0) = std::string (iftext_tmp);
629
630 octave::interpreter& interp = octave::__get_interpreter__ ();
631
632 octave_value_list tmp = interp.feval ("inline", args, 1);
633
634 if (tmp.length () > 0)
635 {
636 retval = tmp(0);
637 return 1;
638 }
639
640#else
641 octave_unused_parameter (loc_id);
642 octave_unused_parameter (name);
643 octave_unused_parameter (retval);
644
645 warn_load ("hdf5");
646#endif
647
648 return -1;
649}
650
651// This function is designed to be passed to H5Giterate, which calls it
652// on each data item in an HDF5 file. For the item whose name is NAME in
653// the group GROUP_ID, this function sets dv->tc to an Octave representation
654// of that item. (dv must be a pointer to hdf5_callback_data.) (It also
655// sets the other fields of dv).
656//
657// It returns 1 on success (in which case H5Giterate stops and returns),
658// -1 on error, and 0 to tell H5Giterate to continue on to the next item
659// (e.g., if NAME was a data type we don't recognize).
660//
661// This function must not throw an exception.
662
663static herr_t
664hdf5_read_next_data_internal (hid_t group_id, const char *name, void *dv)
665{
666 hdf5_callback_data *d = static_cast<hdf5_callback_data *> (dv);
667 hid_t type_id = -1;
668 hid_t type_class_id = -1;
669 hid_t data_id = -1;
670 hid_t subgroup_id = -1;
671 hid_t space_id = -1;
672
673 H5G_stat_t info;
674 herr_t retval = 0;
675 bool ident_valid = octave::valid_identifier (name);
676
677 std::string vname = name;
678
679 octave::type_info& type_info = octave::__get_type_info__ ();
680
681 // Allow identifiers as all digits so we can load lists saved by
682 // earlier versions of Octave.
683
684 if (! ident_valid)
685 {
686 // fix the identifier, replacing invalid chars with underscores
687 vname = make_valid_identifier (vname);
688
689 // check again (in case vname was null, empty, or some such thing):
690 ident_valid = octave::valid_identifier (vname);
691 }
692
693 H5Gget_objinfo (group_id, name, 1, &info);
694
695 if (info.type == H5G_GROUP && ident_valid)
696 {
697#if defined (HAVE_HDF5_18)
698 subgroup_id = H5Gopen (group_id, name, octave_H5P_DEFAULT);
699#else
700 subgroup_id = H5Gopen (group_id, name);
701#endif
702
703 if (subgroup_id < 0)
704 {
705 retval = subgroup_id;
706 goto done;
707 }
708
709 if (hdf5_check_attr (subgroup_id, "OCTAVE_NEW_FORMAT"))
710 {
711#if defined (HAVE_HDF5_18)
712 data_id = H5Dopen (subgroup_id, "type", octave_H5P_DEFAULT);
713#else
714 data_id = H5Dopen (subgroup_id, "type");
715#endif
716
717 if (data_id < 0)
718 {
719 retval = data_id;
720 goto done;
721 }
722
723 type_id = H5Dget_type (data_id);
724
725 type_class_id = H5Tget_class (type_id);
726
727 if (type_class_id != H5T_STRING)
728 goto done;
729
730 space_id = H5Dget_space (data_id);
731 hsize_t rank = H5Sget_simple_extent_ndims (space_id);
732
733 if (rank != 0)
734 goto done;
735
736 int slen = H5Tget_size (type_id);
737 if (slen < 0)
738 goto done;
739
740 OCTAVE_LOCAL_BUFFER (char, typ, slen);
741
742 // create datatype for (null-terminated) string to read into:
743 hid_t st_id = H5Tcopy (H5T_C_S1);
744 H5Tset_size (st_id, slen);
745
746 if (H5Dread (data_id, st_id, octave_H5S_ALL, octave_H5S_ALL,
747 octave_H5P_DEFAULT, typ) < 0)
748 goto done;
749
750 H5Tclose (st_id);
751 H5Dclose (data_id);
752
753 if (std::string (typ, slen-1) == "inline function")
754 {
755 retval = load_inline_fcn (subgroup_id, name, d->tc);
756 }
757 else
758 {
759 d->tc = type_info.lookup_type (std::string (typ, slen-1));
760
761 try
762 {
763 retval = (d->tc.load_hdf5 (subgroup_id, "value") ? 1 : -1);
764 }
765 catch (const octave::execution_exception& ee)
766 {
767 retval = -1;
768 }
769 }
770
771 // check for OCTAVE_GLOBAL attribute:
772 d->global = hdf5_check_attr (subgroup_id, "OCTAVE_GLOBAL");
773
774 H5Gclose (subgroup_id);
775 }
776 else
777 {
778 // It seems that this block only applies to an old list type
779 // and that we shouldn't need to handle the old inline
780 // function type here.
781
782 // an HDF5 group is treated as an octave structure by
783 // default (since that preserves name information), and an
784 // octave list otherwise.
785
786 if (hdf5_check_attr (subgroup_id, "OCTAVE_LIST"))
787 d->tc = type_info.lookup_type ("list");
788 else
789 d->tc = type_info.lookup_type ("struct");
790
791 // check for OCTAVE_GLOBAL attribute:
792 d->global = hdf5_check_attr (subgroup_id, "OCTAVE_GLOBAL");
793
794 H5Gclose (subgroup_id);
795
796 try
797 {
798 retval = (d->tc.load_hdf5 (group_id, name) ? 1 : -1);
799 }
800 catch (const octave::execution_exception& ee)
801 {
802 retval = -1;
803 }
804 }
805
806 }
807 else if (info.type == H5G_DATASET && ident_valid)
808 {
809 // It seems that this block only applies to an old version of
810 // Octave HDF5 files and that it is probably not important to
811 // handle the old inline function type here.
812
813 // For backwards compatibility.
814#if defined (HAVE_HDF5_18)
815 data_id = H5Dopen (group_id, name, octave_H5P_DEFAULT);
816#else
817 data_id = H5Dopen (group_id, name);
818#endif
819
820 if (data_id < 0)
821 {
822 retval = data_id;
823 goto done;
824 }
825
826 type_id = H5Dget_type (data_id);
827
828 type_class_id = H5Tget_class (type_id);
829
830 if (type_class_id == H5T_FLOAT)
831 {
832 space_id = H5Dget_space (data_id);
833
834 hsize_t rank = H5Sget_simple_extent_ndims (space_id);
835
836 if (rank == 0)
837 d->tc = type_info.lookup_type ("scalar");
838 else
839 d->tc = type_info.lookup_type ("matrix");
840
841 H5Sclose (space_id);
842 }
843 else if (type_class_id == H5T_INTEGER)
844 {
845 // What integer type do we really have..
846 std::string int_typ;
847#if defined (HAVE_H5T_GET_NATIVE_TYPE)
848 // FIXME: test this code and activated with an autoconf
849 // test!! It is also incorrect for 64-bit indexing!!
850
851 switch (H5Tget_native_type (type_id, H5T_DIR_ASCEND))
852 {
853 case H5T_NATIVE_CHAR:
854 int_typ = "int8 ";
855 break;
856
857 case H5T_NATIVE_SHORT:
858 int_typ = "int16 ";
859 break;
860
861 case H5T_NATIVE_INT:
862 case H5T_NATIVE_LONG:
863 int_typ = "int32 ";
864 break;
865
866 case H5T_NATIVE_LLONG:
867 int_typ = "int64 ";
868 break;
869
870 case H5T_NATIVE_UCHAR:
871 int_typ = "uint8 ";
872 break;
873
874 case H5T_NATIVE_USHORT:
875 int_typ = "uint16 ";
876 break;
877
878 case H5T_NATIVE_UINT:
879 case H5T_NATIVE_ULONG:
880 int_typ = "uint32 ";
881 break;
882
883 case H5T_NATIVE_ULLONG:
884 int_typ = "uint64 ";
885 break;
886 }
887#else
888 hid_t int_sign = H5Tget_sign (type_id);
889
890 if (int_sign == H5T_SGN_ERROR)
891 warning ("load: can't read '%s' (unknown datatype)", name);
892 else
893 {
894 if (int_sign == H5T_SGN_NONE)
895 int_typ.push_back ('u');
896 int_typ.append ("int");
897
898 int slen = H5Tget_size (type_id);
899 if (slen < 0)
900 warning ("load: can't read '%s' (unknown datatype)", name);
901 else
902 {
903 switch (slen)
904 {
905 case 1:
906 int_typ.append ("8 ");
907 break;
908
909 case 2:
910 int_typ.append ("16 ");
911 break;
912
913 case 4:
914 int_typ.append ("32 ");
915 break;
916
917 case 8:
918 int_typ.append ("64 ");
919 break;
920
921 default:
922 warning ("load: can't read '%s' (unknown datatype)",
923 name);
924 int_typ = "";
925 break;
926 }
927 }
928 }
929#endif
930 if (int_typ == "")
931 warning ("load: can't read '%s' (unknown datatype)", name);
932 else
933 {
934 // Matrix or scalar?
935 space_id = H5Dget_space (data_id);
936
937 hsize_t rank = H5Sget_simple_extent_ndims (space_id);
938
939 if (rank == 0)
940 int_typ.append ("scalar");
941 else
942 int_typ.append ("matrix");
943
944 d->tc = type_info.lookup_type (int_typ);
945 H5Sclose (space_id);
946 }
947 }
948 else if (type_class_id == H5T_STRING)
949 d->tc = type_info.lookup_type ("string");
950 else if (type_class_id == H5T_COMPOUND)
951 {
952 hid_t complex_type = hdf5_make_complex_type (H5T_NATIVE_DOUBLE);
953 hid_t range_type = hdf5_make_range_type (H5T_NATIVE_DOUBLE);
954
955 if (hdf5_types_compatible (type_id, complex_type))
956 {
957 // read complex matrix or scalar variable
958 space_id = H5Dget_space (data_id);
959 hsize_t rank = H5Sget_simple_extent_ndims (space_id);
960
961 if (rank == 0)
962 d->tc = type_info.lookup_type ("complex scalar");
963 else
964 d->tc = type_info.lookup_type ("complex matrix");
965
966 H5Sclose (space_id);
967 }
968 else if (hdf5_types_compatible (type_id, range_type))
969 {
970 // If it's not a complex, check if it's a range
971 d->tc = octave_value_typeinfo::lookup_type ("range");
972 }
973 else // Otherwise, just ignore it with a warning.
974 {
975 warning ("load: can't read '%s' (unknown datatype)", name);
976 retval = 0; // unknown datatype; skip
977 return retval;
978 }
979
980 H5Tclose (range_type);
981 H5Tclose (complex_type);
982 }
983 else
984 {
985 warning ("load: can't read '%s' (unknown datatype)", name);
986 retval = 0; // unknown datatype; skip
987 return retval;
988 }
989
990 // check for OCTAVE_GLOBAL attribute:
991 d->global = hdf5_check_attr (data_id, "OCTAVE_GLOBAL");
992
993 try
994 {
995 retval = (d->tc.load_hdf5 (group_id, name) ? 1 : -1);
996 }
997 catch (const octave::execution_exception& ee)
998 {
999 retval = -1;
1000 }
1001
1002 H5Tclose (type_id);
1003 H5Dclose (data_id);
1004 }
1005
1006 if (! ident_valid)
1007 {
1008 // should we attempt to handle invalid identifiers by converting
1009 // bad characters to '_', say?
1010 warning ("load: skipping invalid identifier '%s' in hdf5 file",
1011 name);
1012 }
1013
1014done:
1015 if (retval < 0)
1016 {
1017 // Must be warning. A call to error aborts and leaves H5Giterate in
1018 // a mangled state that causes segfault on exit (bug #56149).
1019 warning ("load: error while reading hdf5 item '%s'", name);
1020 }
1021
1022 if (retval > 0)
1023 {
1024 // get documentation string, if any:
1025 int comment_length = H5Gget_comment (group_id, name, 0, nullptr);
1026
1027 if (comment_length > 1)
1028 {
1029 OCTAVE_LOCAL_BUFFER (char, tdoc, comment_length);
1030 H5Gget_comment (group_id, name, comment_length, tdoc);
1031 d->doc = tdoc;
1032 }
1033 else if (vname != name)
1034 {
1035 // the name was changed; store the original name
1036 // as the documentation string:
1037 d->doc = name;
1038 }
1039
1040 // copy name (actually, vname):
1041 d->name = vname;
1042 }
1043
1044 return retval;
1045}
1046
1047#endif
1048
1050hdf5_read_next_data (octave_hdf5_id group_id, const char *name, void *dv)
1051{
1052#if defined (HAVE_HDF5)
1053
1054 hid_t new_id = check_hdf5_id_value (group_id, "hdf5_read_next_data");
1055
1056 return hdf5_read_next_data_internal (new_id, name, dv);
1057
1058#else
1059 err_disabled_feature ("hdf5_read_next_data", "HDF5");
1060#endif
1061}
1062
1064hdf5_h5g_iterate (octave_hdf5_id loc_id, const char *name, int *idx,
1065 void *operator_data)
1066{
1067#if defined (HAVE_HDF5)
1068
1069 hid_t new_id = check_hdf5_id_value (loc_id, "hdf5_h5g_iterate");
1070
1071 return H5Giterate (new_id, name, idx, hdf5_read_next_data_internal,
1072 operator_data);
1073
1074#else
1075 err_disabled_feature ("hdf5_h5g_iterate", "HDF5");
1076#endif
1077}
1078
1079// Read the next Octave variable from the stream IS, which must really be an
1080// hdf5_ifstream. Return the variable value in tc, its docstring in doc, and
1081// whether it is global in global. The return value is the name of the
1082// variable, or NULL if none were found or there was an error.
1083std::string
1084read_hdf5_data (std::istream& is, const std::string& /* filename */,
1085 bool& global, octave_value& tc, std::string& doc,
1086 const string_vector& argv, int argv_idx, int argc)
1087{
1088#if defined (HAVE_HDF5)
1089
1090 octave::check_hdf5_types ();
1091
1092 std::string retval;
1093
1094 doc.clear ();
1095
1096 hdf5_ifstream& hs = dynamic_cast<hdf5_ifstream&> (is);
1098
1099 herr_t H5Giterate_retval = -1;
1100
1101 hsize_t num_obj = 0;
1102#if defined (HAVE_HDF5_18)
1103 hid_t group_id = H5Gopen (hs.file_id, "/", octave_H5P_DEFAULT);
1104#else
1105 hid_t group_id = H5Gopen (hs.file_id, "/");
1106#endif
1107 H5Gget_num_objs (group_id, &num_obj);
1108 H5Gclose (group_id);
1109
1110 // For large datasets and out-of-core functionality,
1111 // check if only parts of the data is requested
1112 bool load_named_vars = argv_idx < argc;
1113 while (load_named_vars && hs.current_item < static_cast<int> (num_obj))
1114 {
1115 std::vector<char> var_name;
1116 bool found = false;
1117 std::size_t len = 0;
1118
1119 len = H5Gget_objname_by_idx (hs.file_id, hs.current_item, nullptr, 0);
1120 var_name.resize (len+1);
1121 H5Gget_objname_by_idx (hs.file_id, hs.current_item, &var_name[0], len+1);
1122
1123 for (int i = argv_idx; i < argc; i++)
1124 {
1125 symbol_match pattern (argv[i]);
1126 if (pattern.match (std::string (&var_name[0])))
1127 {
1128 found = true;
1129 break;
1130 }
1131 }
1132
1133 if (found)
1134 break;
1135
1136 hs.current_item++;
1137 }
1138
1139 if (hs.current_item < static_cast<int> (num_obj))
1140 H5Giterate_retval = H5Giterate (hs.file_id, "/", &hs.current_item,
1141 hdf5_read_next_data_internal, &d);
1142
1143 if (H5Giterate_retval > 0)
1144 {
1145 global = d.global;
1146 tc = d.tc;
1147 doc = d.doc;
1148 }
1149 else
1150 {
1151 // An error occurred (H5Giterate_retval < 0),
1152 // or there are no more datasets (H5Giterate_retval == 0).
1153 // hdf5_read_next_data_internal has already printed a warning msg.
1154 }
1155
1156 if (! d.name.empty ())
1157 retval = d.name;
1158
1159 return retval;
1160
1161#else
1162 err_disabled_feature ("read_hdf5_data", "HDF5");
1163#endif
1164}
1165
1166// Add an attribute named attr_name to loc_id (a simple scalar
1167// attribute with value 1). Return value is >= 0 on success.
1169hdf5_add_attr (octave_hdf5_id loc_id, const char *attr_name)
1170{
1171#if defined (HAVE_HDF5)
1172
1173 herr_t retval = 0;
1174
1175 hid_t as_id = H5Screate (H5S_SCALAR);
1176
1177 if (as_id >= 0)
1178 {
1179#if defined (HAVE_HDF5_18)
1180 hid_t a_id = H5Acreate (loc_id, attr_name, H5T_NATIVE_UCHAR,
1182#else
1183 hid_t a_id = H5Acreate (loc_id, attr_name,
1184 H5T_NATIVE_UCHAR, as_id, octave_H5P_DEFAULT);
1185#endif
1186 if (a_id >= 0)
1187 {
1188 unsigned char attr_val = 1;
1189
1190 retval = H5Awrite (a_id, H5T_NATIVE_UCHAR, &attr_val);
1191
1192 H5Aclose (a_id);
1193 }
1194 else
1195 retval = a_id;
1196
1197 H5Sclose (as_id);
1198 }
1199 else
1200 retval = as_id;
1201
1202 return retval;
1203
1204#else
1205 err_disabled_feature ("hdf5_add_attr", "HDF5");
1206#endif
1207}
1208
1211 const char *attr_name, void *buf)
1212{
1213#if defined (HAVE_HDF5)
1214
1215 herr_t retval = 0;
1216
1217 hid_t as_id = H5Screate (H5S_SCALAR);
1218
1219 if (as_id >= 0)
1220 {
1221#if defined (HAVE_HDF5_18)
1222 hid_t a_id = H5Acreate (loc_id, attr_name, type_id,
1224#else
1225 hid_t a_id = H5Acreate (loc_id, attr_name,
1226 type_id, as_id, octave_H5P_DEFAULT);
1227#endif
1228 if (a_id >= 0)
1229 {
1230 retval = H5Awrite (a_id, type_id, buf);
1231
1232 H5Aclose (a_id);
1233 }
1234 else
1235 retval = a_id;
1236
1237 H5Sclose (as_id);
1238 }
1239 else
1240 retval = as_id;
1241
1242 return retval;
1243
1244#else
1245 err_disabled_feature ("hdf5_add_scalar_attr", "HDF5");
1246#endif
1247}
1248
1249// Save an empty matrix, if needed. Returns
1250// > 0 Saved empty matrix
1251// = 0 Not an empty matrix; did nothing
1252// < 0 Error condition
1253int
1254save_hdf5_empty (octave_hdf5_id loc_id, const char *name, const dim_vector& d)
1255{
1256#if defined (HAVE_HDF5)
1257
1258 hsize_t sz = d.length ();
1260 bool empty = false;
1261 hid_t space_hid = -1;
1262 hid_t data_hid = -1;
1263 int retval;
1264 for (hsize_t i = 0; i < sz; i++)
1265 {
1266 dims[i] = d(i);
1267 if (dims[i] < 1)
1268 empty = true;
1269 }
1270
1271 if (! empty)
1272 return 0;
1273
1274 space_hid = H5Screate_simple (1, &sz, nullptr);
1275 if (space_hid < 0) return space_hid;
1276#if defined (HAVE_HDF5_18)
1277 data_hid = H5Dcreate (loc_id, name, H5T_NATIVE_IDX, space_hid,
1280#else
1281 data_hid = H5Dcreate (loc_id, name, H5T_NATIVE_IDX, space_hid,
1283#endif
1284 if (data_hid < 0)
1285 {
1286 H5Sclose (space_hid);
1287 return data_hid;
1288 }
1289
1290 retval = H5Dwrite (data_hid, H5T_NATIVE_IDX, octave_H5S_ALL, octave_H5S_ALL,
1291 octave_H5P_DEFAULT, dims) >= 0;
1292
1293 H5Dclose (data_hid);
1294 H5Sclose (space_hid);
1295
1296 if (retval)
1297 retval = hdf5_add_attr (loc_id, "OCTAVE_EMPTY_MATRIX");
1298
1299 return (retval == 0 ? 1 : retval);
1300
1301#else
1302 err_disabled_feature ("save_hdf5_empty", "HDF5");
1303#endif
1304}
1305
1306// Load an empty matrix, if needed. Returns
1307// > 0 loaded empty matrix, dimensions returned
1308// = 0 Not an empty matrix; did nothing
1309// < 0 Error condition
1310int
1311load_hdf5_empty (octave_hdf5_id loc_id, const char *name, dim_vector& d)
1312{
1313#if defined (HAVE_HDF5)
1314
1315 if (! hdf5_check_attr (loc_id, "OCTAVE_EMPTY_MATRIX"))
1316 return 0;
1317
1318 hsize_t hdims, maxdims;
1319#if defined (HAVE_HDF5_18)
1320 hid_t data_hid = H5Dopen (loc_id, name, octave_H5P_DEFAULT);
1321#else
1322 hid_t data_hid = H5Dopen (loc_id, name);
1323#endif
1324 hid_t space_id = H5Dget_space (data_hid);
1325 H5Sget_simple_extent_dims (space_id, &hdims, &maxdims);
1326 int retval;
1327
1328 OCTAVE_LOCAL_BUFFER (octave_idx_type, dims, hdims);
1329
1330 retval = H5Dread (data_hid, H5T_NATIVE_IDX, octave_H5S_ALL, octave_H5S_ALL,
1331 octave_H5P_DEFAULT, dims);
1332 if (retval >= 0)
1333 {
1334 d.resize (hdims);
1335 for (hsize_t i = 0; i < hdims; i++)
1336 d(i) = dims[i];
1337 }
1338
1339 H5Sclose (space_id);
1340 H5Dclose (data_hid);
1341
1342 return (retval == 0 ? hdims : retval);
1343
1344#else
1345 err_disabled_feature ("load_hdf5_empty", "HDF5");
1346#endif
1347}
1348
1349// save_type_to_hdf5 is not currently used, since hdf5 doesn't yet support
1350// automatic float<->integer conversions:
1351
1352// return the HDF5 type id corresponding to the Octave save_type
1353
1356{
1357#if defined (HAVE_HDF5)
1358# if defined (HAVE_HDF5_INT2FLOAT_CONVERSIONS)
1359
1360 switch (st)
1361 {
1362 case LS_U_CHAR:
1363 return H5T_NATIVE_UCHAR;
1364
1365 case LS_U_SHORT:
1366 return H5T_NATIVE_USHORT;
1367
1368 case LS_U_INT:
1369 return H5T_NATIVE_UINT;
1370
1371 case LS_CHAR:
1372 return H5T_NATIVE_CHAR;
1373
1374 case LS_SHORT:
1375 return H5T_NATIVE_SHORT;
1376
1377 case LS_INT:
1378 return H5T_NATIVE_INT;
1379
1380 case LS_FLOAT:
1381 return H5T_NATIVE_FLOAT;
1382
1383 case LS_DOUBLE:
1384 default:
1385 return H5T_NATIVE_DOUBLE;
1386 }
1387
1388# else
1389
1390 octave_unused_parameter (st);
1391
1392 return -1;
1393
1394# endif
1395
1396#else
1397
1398 octave_unused_parameter (st);
1399
1400 err_disabled_feature ("save_type_to_hdf5", "HDF5");
1401
1402#endif
1403}
1404
1405// Add the data from TC to the HDF5 location loc_id, which could
1406// be either a file or a group within a file. Return true if
1407// successful. This function calls itself recursively for lists
1408// (stored as HDF5 groups).
1409
1410bool
1412 const std::string& name, const std::string& doc,
1413 bool mark_global, bool save_as_floats)
1414{
1415#if defined (HAVE_HDF5)
1416
1417 hsize_t dims[3];
1418 hid_t type_id, space_id, data_id, data_type_id;
1419 type_id = space_id = data_id = data_type_id = -1;
1420
1421 bool retval = false;
1422 octave_value val = tc;
1423 // FIXME: diagonal & permutation matrices currently don't know how to save
1424 // themselves, so we convert them first to normal matrices using A = A(:,:).
1425 // This is a temporary hack.
1426 if (val.is_diag_matrix () || val.is_perm_matrix ()
1428 val = val.full_value ();
1429
1430 std::string t = val.type_name ();
1431#if defined (HAVE_HDF5_18)
1432 data_id = H5Gcreate (loc_id, name.c_str (), octave_H5P_DEFAULT,
1434#else
1435 data_id = H5Gcreate (loc_id, name.c_str (), 0);
1436#endif
1437 if (data_id < 0)
1438 goto error_cleanup;
1439
1440 // attach the type of the variable
1441 type_id = H5Tcopy (H5T_C_S1);
1442 H5Tset_size (type_id, t.length () + 1);
1443 if (type_id < 0)
1444 goto error_cleanup;
1445
1446 dims[0] = 0;
1447 space_id = H5Screate_simple (0, dims, nullptr);
1448 if (space_id < 0)
1449 goto error_cleanup;
1450#if defined (HAVE_HDF5_18)
1451 data_type_id = H5Dcreate (data_id, "type", type_id, space_id,
1454#else
1455 data_type_id = H5Dcreate (data_id, "type", type_id, space_id,
1457#endif
1458 if (data_type_id < 0
1459 || H5Dwrite (data_type_id, type_id, octave_H5S_ALL, octave_H5S_ALL,
1460 octave_H5P_DEFAULT, t.c_str ()) < 0)
1461 goto error_cleanup;
1462
1463 // Now call the real function to save the variable
1464 retval = val.save_hdf5 (data_id, "value", save_as_floats);
1465
1466 // attach doc string as comment:
1467 if (retval && doc.length () > 0
1468 && H5Gset_comment (loc_id, name.c_str (), doc.c_str ()) < 0)
1469 retval = false;
1470
1471 // if it's global, add an attribute "OCTAVE_GLOBAL" with value 1
1472 if (retval && mark_global)
1473 retval = hdf5_add_attr (data_id, "OCTAVE_GLOBAL") >= 0;
1474
1475 // We are saving in the new variable format, so mark it
1476 if (retval)
1477 retval = hdf5_add_attr (data_id, "OCTAVE_NEW_FORMAT") >= 0;
1478
1479error_cleanup:
1480
1481 if (data_type_id >= 0)
1482 H5Dclose (data_type_id);
1483
1484 if (type_id >= 0)
1485 H5Tclose (type_id);
1486
1487 if (space_id >= 0)
1488 H5Sclose (space_id);
1489
1490 if (data_id >= 0)
1491 H5Gclose (data_id);
1492
1493 if (! retval)
1494 error ("save: error while writing '%s' to hdf5 file", name.c_str ());
1495
1496 return retval;
1497
1498#else
1499 err_disabled_feature ("add_hdf5_data", "HDF5");
1500#endif
1501}
1502
1503// Write data from TC in HDF5 (binary) format to the stream OS,
1504// which must be an hdf5_ofstream, returning true on success.
1505
1506bool
1507save_hdf5_data (std::ostream& os, const octave_value& tc,
1508 const std::string& name, const std::string& doc,
1509 bool mark_global, bool save_as_floats)
1510{
1511#if defined (HAVE_HDF5)
1512
1513 octave::check_hdf5_types ();
1514
1515 hdf5_ofstream& hs = dynamic_cast<hdf5_ofstream&> (os);
1516
1517 return add_hdf5_data (hs.file_id, tc, name, doc,
1518 mark_global, save_as_floats);
1519
1520#else
1521 err_disabled_feature ("save_hdf5_data", "HDF5");
1522#endif
1523}
1524
1525#endif
Vector representing the dimensions (size) of an Array.
Definition dim-vector.h:90
void open(const char *name, int mode, int)
Definition ls-hdf5.cc:130
void open_create(const char *name, int mode)
Definition ls-hdf5.cc:148
octave_hdf5_id file_id
Definition ls-hdf5.h:47
static int static_type_id()
octave_idx_type length() const
Definition ovl.h:111
bool is_diag_matrix() const
Definition ov.h:631
bool is_perm_matrix() const
Definition ov.h:634
bool save_hdf5(octave_hdf5_id loc_id, const char *name, bool save_as_floats)
Definition ov.h:1391
int type_id() const
Definition ov.h:1358
octave_value full_value() const
Definition ov.h:435
std::string type_name() const
Definition ov.h:1360
bool match(const std::string &sym)
Definition glob-match.cc:86
octave_value lookup_type(const std::string &nm)
const octave_hdf5_id octave_H5P_DEFAULT
const octave_hdf5_id octave_H5E_DEFAULT
const octave_hdf5_id octave_H5S_ALL
save_type
Definition data-conv.h:85
@ LS_U_CHAR
Definition data-conv.h:86
@ LS_DOUBLE
Definition data-conv.h:93
@ LS_U_SHORT
Definition data-conv.h:87
@ LS_FLOAT
Definition data-conv.h:92
@ LS_SHORT
Definition data-conv.h:90
@ LS_CHAR
Definition data-conv.h:89
@ LS_INT
Definition data-conv.h:91
@ LS_U_INT
Definition data-conv.h:88
void warning(const char *fmt,...)
Definition error.cc:1078
void error(const char *fmt,...)
Definition error.cc:1003
void err_disabled_feature(const std::string &fcn, const std::string &feature, const std::string &pkg)
Definition errwarn.cc:53
F77_RET_T const F77_DBLE const F77_DBLE F77_DBLE * d
bool save_hdf5_data(std::ostream &os, const octave_value &tc, const std::string &name, const std::string &doc, bool mark_global, bool save_as_floats)
Definition ls-hdf5.cc:1507
int save_hdf5_empty(octave_hdf5_id loc_id, const char *name, const dim_vector &d)
Definition ls-hdf5.cc:1254
std::string read_hdf5_data(std::istream &is, const std::string &, bool &global, octave_value &tc, std::string &doc, const string_vector &argv, int argv_idx, int argc)
Definition ls-hdf5.cc:1084
octave_hdf5_err hdf5_add_attr(octave_hdf5_id loc_id, const char *attr_name)
Definition ls-hdf5.cc:1169
bool hdf5_get_scalar_attr(octave_hdf5_id loc_id, octave_hdf5_id type_id, const char *attr_name, void *buf)
Definition ls-hdf5.cc:352
int load_hdf5_empty(octave_hdf5_id loc_id, const char *name, dim_vector &d)
Definition ls-hdf5.cc:1311
octave_hdf5_id save_type_to_hdf5(save_type st)
Definition ls-hdf5.cc:1355
octave_hdf5_id hdf5_make_complex_type(octave_hdf5_id num_type)
Definition ls-hdf5.cc:410
octave_hdf5_err hdf5_read_next_data(octave_hdf5_id group_id, const char *name, void *dv)
Definition ls-hdf5.cc:1050
octave_hdf5_err hdf5_h5g_iterate(octave_hdf5_id loc_id, const char *name, int *idx, void *operator_data)
Definition ls-hdf5.cc:1064
bool hdf5_types_compatible(octave_hdf5_id t1, octave_hdf5_id t2)
Definition ls-hdf5.cc:274
octave_hdf5_err hdf5_add_scalar_attr(octave_hdf5_id loc_id, octave_hdf5_id type_id, const char *attr_name, void *buf)
Definition ls-hdf5.cc:1210
bool add_hdf5_data(octave_hdf5_id loc_id, const octave_value &tc, const std::string &name, const std::string &doc, bool mark_global, bool save_as_floats)
Definition ls-hdf5.cc:1411
bool hdf5_check_attr(octave_hdf5_id loc_id, const char *attr_name)
Definition ls-hdf5.cc:305
int octave_hdf5_err
int64_t octave_hdf5_id
#define H5T_NATIVE_IDX
Definition oct-hdf5.h:42
#define OCTAVE_LOCAL_BUFFER(T, buf, size)
Definition oct-locbuf.h:44
F77_RET_T len
Definition xerbla.cc:61