bugGNU Octave - Bugs: bug #58268, hdf5 load crashes on strings

 
 

bug #58268: hdf5 load crashes on strings

Submitter:  Christian Häggström <saturn>
Submitted:  Tue 28 Apr 2020 09:46:57 PM UTC
   
 
Category:  Octave Function Severity:  3 - Normal
Priority:  5 - Normal Item Group:  Segfault, Bus Error, etc.
Status:  Fixed Assigned to:  None
Originator Name:  Open/Closed:  * Closed
Release:  * 5.2.0 Operating System:  * GNU/Linux
Fixed Release:  None Planned Release:  None
* Mandatory Fields

Add a New Comment Rich Markup
   

Jump to the original submission

Wed 29 Apr 2020 09:46:06 AM UTC, comment #7: 

Can confirm the fix works as advertised. The valgrind output is clean.

(My previous mistake was forgetting to run hg update tip after hg pull)

Elvis Stansvik <estan>
Wed 29 Apr 2020 09:34:55 AM UTC, comment #6: 

Sorry, my bad! Somehow I did not have your patch applied after all (I'm Mercurial newbie). I will correct the mistake and report back.

Elvis Stansvik <estan>
Wed 29 Apr 2020 09:25:54 AM UTC, comment #5: 

Valgrind reports:


$ valgrind octave-cli --eval "load outfile.hdf5"
==24612== Memcheck, a memory error detector
==24612== Copyright (C) 2002-2017, and GNU GPL'd, by Julian Seward et al.
==24612== Using Valgrind-3.15.0 and LibVEX; rerun with -h for copyright info
==24612== Command: octave-cli --eval load\ outfile.hdf5
==24612==
==24612== Invalid write of size 1
==24612==    at 0x4840273: memcpy@GLIBC_2.2.5 (in /usr/lib/x86_64-linux-gnu/valgrind/vgpreload_memcheck-amd64-linux.so)
==24612==    by 0x737A37F: memcpy (string_fortified.h:34)
==24612==    by 0x737A37F: H5D__scatter_mem (H5Dscatgath.c:340)
==24612==    by 0x737A37F: H5D__scatter_mem (H5Dscatgath.c:291)
==24612==    by 0x737A9C9: H5D__scatgath_read (H5Dscatgath.c:567)
==24612==    by 0x7362FA3: H5D__contig_read (H5Dcontig.c:595)
==24612==    by 0x7376513: H5D__read (H5Dio.c:600)
==24612==    by 0x7376A0C: H5Dread (H5Dio.c:198)
==24612==    by 0x553EE84: octave_char_matrix_str::load_hdf5(long, char const*) (ov-str-mat.cc:703)
==24612==    by 0x5A52675: octave_value::load_hdf5(long, char const*) (ov.h:1321)
==24612==    by 0x5A51008: hdf5_read_next_data_internal(long, char const*, void*) (ls-hdf5.cc:748)
==24612==    by 0x73E7835: H5G_iterate_cb (H5Gint.c:818)
==24612==    by 0x73EF237: H5G__node_iterate (H5Gnode.c:1002)
==24612==    by 0x7320FC5: H5B__iterate_helper (H5B.c:1165)
==24612==  Address 0xfbf7ae8 is 0 bytes after a block of size 24 alloc'd
==24612==    at 0x483C583: operator new[](unsigned long) (in /usr/lib/x86_64-linux-gnu/valgrind/vgpreload_memcheck-amd64-linux.so)
==24612==    by 0x546FB2C: std::_MakeUniq<char []>::__array std::make_unique<char []>(unsigned long) (unique_ptr.h:863)
==24612==    by 0x553EE00: octave_char_matrix_str::load_hdf5(long, char const*) (ov-str-mat.cc:698)
==24612==    by 0x5A52675: octave_value::load_hdf5(long, char const*) (ov.h:1321)
==24612==    by 0x5A51008: hdf5_read_next_data_internal(long, char const*, void*) (ls-hdf5.cc:748)
==24612==    by 0x73E7835: H5G_iterate_cb (H5Gint.c:818)
==24612==    by 0x73EF237: H5G__node_iterate (H5Gnode.c:1002)
==24612==    by 0x7320FC5: H5B__iterate_helper (H5B.c:1165)
==24612==    by 0x732255D: H5B_iterate (H5B.c:1210)
==24612==    by 0x73F546C: H5G__stab_iterate (H5Gstab.c:556)
==24612==    by 0x73F1E68: H5G__obj_iterate (H5Gobj.c:697)
==24612==    by 0x73E8BCE: H5G_iterate (H5Gint.c:892)
==24612==
==24612== Invalid read of size 1
==24612==    at 0x483EFB4: __strlen_sse2 (in /usr/lib/x86_64-linux-gnu/valgrind/vgpreload_memcheck-amd64-linux.so)
==24612==    by 0x66626B1: charNDArray::charNDArray(char const*) (chNDArray.cc:55)
==24612==    by 0x553F6C8: charMatrix::charMatrix(char const*) (chMatrix.h:70)
==24612==    by 0x553EEDE: octave_char_matrix_str::load_hdf5(long, char const*) (ov-str-mat.cc:713)
==24612==    by 0x5A52675: octave_value::load_hdf5(long, char const*) (ov.h:1321)
==24612==    by 0x5A51008: hdf5_read_next_data_internal(long, char const*, void*) (ls-hdf5.cc:748)
==24612==    by 0x73E7835: H5G_iterate_cb (H5Gint.c:818)
==24612==    by 0x73EF237: H5G__node_iterate (H5Gnode.c:1002)
==24612==    by 0x7320FC5: H5B__iterate_helper (H5B.c:1165)
==24612==    by 0x732255D: H5B_iterate (H5B.c:1210)
==24612==    by 0x73F546C: H5G__stab_iterate (H5Gstab.c:556)
==24612==    by 0x73F1E68: H5G__obj_iterate (H5Gobj.c:697)
==24612==  Address 0xfbf7ae8 is 0 bytes after a block of size 24 alloc'd
==24612==    at 0x483C583: operator new[](unsigned long) (in /usr/lib/x86_64-linux-gnu/valgrind/vgpreload_memcheck-amd64-linux.so)
==24612==    by 0x546FB2C: std::_MakeUniq<char []>::__array std::make_unique<char []>(unsigned long) (unique_ptr.h:863)
==24612==    by 0x553EE00: octave_char_matrix_str::load_hdf5(long, char const*) (ov-str-mat.cc:698)
==24612==    by 0x5A52675: octave_value::load_hdf5(long, char const*) (ov.h:1321)
==24612==    by 0x5A51008: hdf5_read_next_data_internal(long, char const*, void*) (ls-hdf5.cc:748)
==24612==    by 0x73E7835: H5G_iterate_cb (H5Gint.c:818)
==24612==    by 0x73EF237: H5G__node_iterate (H5Gnode.c:1002)
==24612==    by 0x7320FC5: H5B__iterate_helper (H5B.c:1165)
==24612==    by 0x732255D: H5B_iterate (H5B.c:1210)
==24612==    by 0x73F546C: H5G__stab_iterate (H5Gstab.c:556)
==24612==    by 0x73F1E68: H5G__obj_iterate (H5Gobj.c:697)
==24612==    by 0x73E8BCE: H5G_iterate (H5Gint.c:892)
==24612==
==24612==
==24612== HEAP SUMMARY:
==24612==     in use at exit: 261,635 bytes in 4,145 blocks
==24612==   total heap usage: 186,092 allocs, 181,947 frees, 31,890,069 bytes allocated
==24612==
==24612== LEAK SUMMARY:
==24612==    definitely lost: 935 bytes in 61 blocks
==24612==    indirectly lost: 159,921 bytes in 2,727 blocks
==24612==      possibly lost: 688 bytes in 9 blocks
==24612==    still reachable: 100,091 bytes in 1,348 blocks
==24612==                       of which reachable via heuristic:
==24612==                         newarray           : 424 bytes in 1 blocks
==24612==         suppressed: 0 bytes in 0 blocks
==24612== Rerun with --leak-check=full to see details of leaked memory
==24612==
==24612== For lists of detected and suppressed errors, rerun with: -s
==24612== ERROR SUMMARY: 2 errors from 2 contexts (suppressed: 0 from 0)
$


Elvis Stansvik <estan>
Wed 29 Apr 2020 09:13:39 AM UTC, comment #4: 

Thanks for the fix Rik! However, I tried the test case out, using Octave from hg and with debugging symbols, and I still get a crash, though with the following backtrace:


(gdb) run --eval "load outfile.hdf5"
Starting program: /home/buildbot/my_octave/bin/octave-cli --eval "load outfile.hdf5"
[Thread debugging using libthread_db enabled]
Using host libthread_db library "/lib/x86_64-linux-gnu/libthread_db.so.1".
[New Thread 0x7fffef6ba700 (LWP 24551)]
free(): invalid pointer

Thread 1 "octave-cli" received signal SIGABRT, Aborted.
__GI_raise (sig=sig@entry=6) at ../sysdeps/unix/sysv/linux/raise.c:50
50      ../sysdeps/unix/sysv/linux/raise.c: No such file or directory.
(gdb) bt
#0  __GI_raise (sig=sig@entry=6) at ../sysdeps/unix/sysv/linux/raise.c:50
#1  0x00007ffff5574859 in __GI_abort () at abort.c:79
#2  0x00007ffff55df3ee in __libc_message (action=action@entry=do_abort, fmt=fmt@entry=0x7ffff5709285 "%s\n") at ../sysdeps/posix/libc_fatal.c:155
#3  0x00007ffff55e747c in malloc_printerr (str=str@entry=0x7ffff57074ae "free(): invalid pointer") at malloc.c:5347
#4  0x00007ffff55e8cac in _int_free (av=<optimized out>, p=<optimized out>, have_lock=0) at malloc.c:4173
#5  0x00007ffff7646671 in octave::tree_constant::~tree_constant (this=0x55555597e530, __in_chrg=<optimized out>) at ./../libinterp/parse-tree/pt-const.h:70
#6  0x00007ffff763bb5b in octave::tree_simple_assignment::~tree_simple_assignment (this=0x55555597e430, __in_chrg=<optimized out>) at ./../libinterp/parse-tree/pt-assign.cc:56
#7  0x00007ffff763bb86 in octave::tree_simple_assignment::~tree_simple_assignment (this=0x55555597e430, __in_chrg=<optimized out>) at ./../libinterp/parse-tree/pt-assign.cc:57
#8  0x00007ffff767bd76 in octave::tree_statement::~tree_statement (this=0x55555597e5d0, __in_chrg=<optimized out>) at ./../libinterp/parse-tree/pt-stmt.cc:59
#9  0x00007ffff767bdbe in octave::tree_statement::~tree_statement (this=0x55555597e5d0, __in_chrg=<optimized out>) at ./../libinterp/parse-tree/pt-stmt.cc:61
#10 0x00007ffff7623dc3 in octave::tree_statement_list::~tree_statement_list (this=0x555555859940, __in_chrg=<optimized out>) at ./../libinterp/parse-tree/pt-stmt.h:166
#11 0x00007ffff7623e18 in octave::tree_statement_list::~tree_statement_list (this=0x555555859940, __in_chrg=<optimized out>) at ./../libinterp/parse-tree/pt-stmt.h:169
#12 0x00007ffff75577af in octave_user_code::~octave_user_code (this=0x55555599c350, __in_chrg=<optimized out>) at ./../libinterp/octave-value/ov-usr-fcn.cc:84
#13 0x00007ffff7558b86 in octave_user_function::~octave_user_function (this=0x55555599c350, __in_chrg=<optimized out>) at ./../libinterp/octave-value/ov-usr-fcn.cc:220
#14 0x00007ffff7558ba6 in octave_user_function::~octave_user_function (this=0x55555599c350, __in_chrg=<optimized out>) at ./../libinterp/octave-value/ov-usr-fcn.cc:230
#15 0x00007ffff7043951 in octave_value::operator= (this=0x5555557f52e0, a=...) at ./../libinterp/octave-value/ov.h:393
#16 0x00007ffff7bad9cd in octave::fcn_info::fcn_info_rep::clear_user_function (this=0x5555557f51d0, force=false) at ./../libinterp/corefcn/fcn-info.h:162
#17 0x00007ffff7badb9f in octave::fcn_info::fcn_info_rep::clear (this=0x5555557f51d0, force=false) at ./../libinterp/corefcn/fcn-info.h:187
#18 0x00007ffff7badfe8 in octave::fcn_info::clear (this=0x5555559433a0, force=false) at ./../libinterp/corefcn/fcn-info.h:322
#19 0x00007ffff7baaf06 in octave::symbol_table::clear_functions (this=0x5555555d0198, force=false) at ./../libinterp/corefcn/symtab.cc:339
#20 0x00007ffff7bac0cd in octave::symbol_table::cleanup (this=0x5555555d0198) at ./../libinterp/corefcn/symtab.cc:603
#21 0x00007ffff7a180f7 in octave::interpreter::cleanup (this=0x5555555cf4a0) at ./../libinterp/corefcn/interpreter.cc:1269
#22 0x00007ffff7a155c1 in octave::interpreter::~interpreter (this=0x5555555cf4a0, __in_chrg=<optimized out>) at ./../libinterp/corefcn/interpreter.cc:638
#23 0x00007ffff702946f in octave::application::delete_interpreter (this=0x7fffffffe170) at ./../libinterp/octave.cc:320
#24 0x00007ffff7029841 in octave::cli_application::execute (this=0x7fffffffe170) at ./../libinterp/octave.cc:381
#25 0x00005555555566ec in main (argc=3, argv=0x7fffffffe478) at ./../src/main-cli.cc:95
(gdb)


This was with changeset be19ade4acc3:


(env) buildbot@759e20843d04:~/src/octave$ hg id -i
be19ade4acc3
(env) buildbot@759e20843d04:~/src/octave$


Elvis Stansvik <estan>
Wed 29 Apr 2020 05:54:15 AM UTC, comment #3: 

Thanks Rik for a swift assessment and a fix! Much appreciated.

Christian Häggström <saturn>
Tue 28 Apr 2020 11:26:10 PM UTC, comment #2: 

The problem was that in one instance, the output buffer was sized to strlen rather than strlen+1 to accommodate the null byte termination.  I fixed that in this changeset https://hg.savannah.gnu.org/hgweb/octave/rev/9646d752c76c.  This will be a part of the next upcoming release of Octave which is version 6.1.

Rik <rik5>
Group administrator
Tue 28 Apr 2020 10:39:01 PM UTC, comment #1: 

Because Octave did not create the file it is unlikely to load correctly.  But, as you point out this shouldn't, cause an error.  Working with the development branch of Octave where I have debugging symbols the file loads.  However, the very next command, whatever it is, causes a segfault.  The backtrace doesn't show anything related to hdf5 which I think indicates that the H5Dread routine itself is causing memory corruption.

Rik <rik5>
Group administrator
Tue 28 Apr 2020 09:46:57 PM UTC, original submission:  

Minimal testcase can be created with this Python snippet, also attached.

python3 -c "import h5py, numpy; f = h5py.File('outfile.hdf5', 'w'); f.create_dataset('scan/date', data=numpy.string_('Mon Jun 18 03:34:30 2018')); f.close()"



$ octave-cli --version
GNU Octave, version 5.2.0

$ octave-cli --eval "load outfile.hdf5"
free(): invalid next size (fast)
fatal: caught signal Aborted -- stopping myself...
Aborted

$ valgrind octave-cli --eval "load outfile.hdf5"
==238738== Invalid write of size 1
==238738==    at 0x483B114: memcpy@GLIBC_2.2.5 (vg_replace_strmem.c:1034)
==238738==    by 0x715A2CF: H5D__scatter_mem (in /usr/lib/x86_64-linux-gnu/libhdf5_serial.so.103.2.0)
==238738==    by 0x715A911: H5D__scatgath_read (in /usr/lib/x86_64-linux-gnu/libhdf5_serial.so.103.2.0)
==238738==    by 0x7142E22: H5D__contig_read (in /usr/lib/x86_64-linux-gnu/libhdf5_serial.so.103.2.0)
==238738==    by 0x7156755: H5D__read (in /usr/lib/x86_64-linux-gnu/libhdf5_serial.so.103.2.0)
==238738==    by 0x7156B68: H5Dread (in /usr/lib/x86_64-linux-gnu/libhdf5_serial.so.103.2.0)
==238738==    by 0x5217072: octave_char_matrix_str::load_hdf5(long, char const*) (in /usr/lib/x86_64-linux-gnu/liboctinterp.so.7.0.1)
==238738==    by 0x56ACBBB: ??? (in /usr/lib/x86_64-linux-gnu/liboctinterp.so.7.0.1)
==238738==    by 0x71C44ED: ??? (in /usr/lib/x86_64-linux-gnu/libhdf5_serial.so.103.2.0)
==238738==    by 0x71CBAC7: H5G__node_iterate (in /usr/lib/x86_64-linux-gnu/libhdf5_serial.so.103.2.0)
==238738==    by 0x70FEA25: ??? (in /usr/lib/x86_64-linux-gnu/libhdf5_serial.so.103.2.0)
==238738==    by 0x70FFFA9: H5B_iterate (in /usr/lib/x86_64-linux-gnu/libhdf5_serial.so.103.2.0)
==238738==  Address 0xd124428 is 0 bytes after a block of size 24 alloc'd
==238738==    at 0x483750F: operator new[](unsigned long) (vg_replace_malloc.c:433)
==238738==    by 0x521700F: octave_char_matrix_str::load_hdf5(long, char const*) (in /usr/lib/x86_64-linux-gnu/liboctinterp.so.7.0.1)

It looks to me that octave forgets to allocate space for the trailing NUL byte. Please observe that this file was not created by Octave itself, but even if Octave rejects the file it should not crash.

Christian Häggström <saturn>

 

(Note: upload size limit is set to 16384 kB, after insertion of the required escape characters.)

Attach Files:
   
   
Comment:
   

Attached Files
file #48962:  outfile.hdf5 added by saturn (5KiB - application/octet-stream)

 

Depends on the following items: None found

Items that depend on this one: None found

 

Carbon-Copy List
  • -email is unavailable- added by estan (Posted a comment)
  • -email is unavailable- added by rik5 (Posted a comment)
  • -email is unavailable- added by saturn (Submitted the item)
  •  

    There are 0 votes so far. Votes easily highlight which items people would like to see resolved in priority, independently of the priority of the item set by tracker managers.

    Only group members can vote.

     

    Follow 4 latest changes.

    Date Changed by Updated Field Previous Value => Replaced by
    2020-04-28 rik5 StatusConfirmed Fixed
        Open/ClosedOpen Closed
    2020-04-28 rik5 StatusNone Confirmed
    2020-04-28 saturn Attached File- Added outfile.hdf5, #48962

    Back to the top

    Powered by Savane 3.13-caa5.
    Corresponding source code