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