summaryrefslogtreecommitdiffstats
path: root/src/armadillo/include/armadillo_bits/hdf5_misc.hpp
blob: 0dd4a7a18882fb46ca0559b75d3818e644dd3922 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
// SPDX-License-Identifier: Apache-2.0
// 
// Copyright 2008-2016 Conrad Sanderson (http://conradsanderson.id.au)
// Copyright 2008-2016 National ICT Australia (NICTA)
// 
// Licensed under the Apache License, Version 2.0 (the "License");
// you may not use this file except in compliance with the License.
// You may obtain a copy of the License at
// http://www.apache.org/licenses/LICENSE-2.0
// 
// Unless required by applicable law or agreed to in writing, software
// distributed under the License is distributed on an "AS IS" BASIS,
// WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
// See the License for the specific language governing permissions and
// limitations under the License.
// ------------------------------------------------------------------------


//! \addtogroup hdf5_misc
//! @{


#if defined(ARMA_USE_HDF5)
namespace hdf5_misc
{


//! Given a certain type, find the corresponding HDF5 datatype.  This can't be
//! done entirely at compile time, unfortunately, because the H5T_* macros
//! depend on function calls.
template< typename eT >
inline
hid_t
get_hdf5_type()
  {
  return -1; // Return invalid.
  }



//! Specializations for each valid element type
//! (taken from all the possible typedefs of {u8, s8, ..., u64, s64} and the other native types.  
//! We can't use the actual u8/s8 typedefs because their relations to the H5T_... types are unclear.
template<>
inline
hid_t
get_hdf5_type< unsigned char >()
  {
  return H5Tcopy(H5T_NATIVE_UCHAR);
  }

template<>
inline
hid_t
get_hdf5_type< char >()
  {
  return H5Tcopy(H5T_NATIVE_CHAR);
  }

template<>
inline
hid_t
get_hdf5_type< short >()
  {
  return H5Tcopy(H5T_NATIVE_SHORT);
  }

template<>
inline
hid_t
get_hdf5_type< unsigned short >()
  {
  return H5Tcopy(H5T_NATIVE_USHORT);
  }

template<>
inline
hid_t
get_hdf5_type< int >()
  {
  return H5Tcopy(H5T_NATIVE_INT);
  }

template<>
inline
hid_t
get_hdf5_type< unsigned int >()
  {
  return H5Tcopy(H5T_NATIVE_UINT);
  }

template<>
inline
hid_t
get_hdf5_type< long >()
  {
  return H5Tcopy(H5T_NATIVE_LONG);
  }

template<>
inline
hid_t
get_hdf5_type< unsigned long >()
  {
  return H5Tcopy(H5T_NATIVE_ULONG);
  }

template<>
inline
hid_t
get_hdf5_type< long long >()
  {
  return H5Tcopy(H5T_NATIVE_LLONG);
  }

template<>
inline
hid_t
get_hdf5_type< unsigned long long >()
  {
  return H5Tcopy(H5T_NATIVE_ULLONG);
  }

template<>
inline
hid_t
get_hdf5_type< float >()
  {
  return H5Tcopy(H5T_NATIVE_FLOAT);
  }

template<>
inline
hid_t
get_hdf5_type< double >()
  {
  return H5Tcopy(H5T_NATIVE_DOUBLE);
  }



//! Utility hid_t since HOFFSET() won't work with std::complex.
template<typename eT>
struct hdf5_complex_t
  {
  eT real;
  eT imag;
  };



template<>
inline
hid_t
get_hdf5_type< std::complex<float> >()
  {
  hid_t type = H5Tcreate(H5T_COMPOUND, sizeof(hdf5_complex_t<float>));
  
  H5Tinsert(type, "real", HOFFSET(hdf5_complex_t<float>, real), H5T_NATIVE_FLOAT);
  H5Tinsert(type, "imag", HOFFSET(hdf5_complex_t<float>, imag), H5T_NATIVE_FLOAT);
  
  return type;
  }



template<>
inline
hid_t
get_hdf5_type< std::complex<double> >()
  {
  hid_t type = H5Tcreate(H5T_COMPOUND, sizeof(hdf5_complex_t<double>));

  H5Tinsert(type, "real", HOFFSET(hdf5_complex_t<double>, real), H5T_NATIVE_DOUBLE);
  H5Tinsert(type, "imag", HOFFSET(hdf5_complex_t<double>, imag), H5T_NATIVE_DOUBLE);

  return type;
  }



// Compare datatype against all supported types.
inline
bool
is_supported_arma_hdf5_type(hid_t datatype)
  {
  hid_t  search_type;
  
  bool is_equal;
  
  
  // start with most likely used types: double, complex<double>, float, complex<float>
  
  search_type = get_hdf5_type<double>();
  is_equal = ( H5Tequal(datatype, search_type) > 0 );
  H5Tclose(search_type);
  if(is_equal) { return true; }
  
  search_type = get_hdf5_type< std::complex<double> >();
  is_equal = ( H5Tequal(datatype, search_type) > 0 );
  H5Tclose(search_type);
  if(is_equal) { return true; }
  
  search_type = get_hdf5_type<float>();
  is_equal = ( H5Tequal(datatype, search_type) > 0 );
  H5Tclose(search_type);
  if(is_equal) { return true; }
  
  search_type = get_hdf5_type< std::complex<float> >();
  is_equal = ( H5Tequal(datatype, search_type) > 0 );
  H5Tclose(search_type);
  if(is_equal) { return true; }
  
  
  // remaining supported types: u8, s8, u16, s16, u32, s32, u64, s64, ulng_t, slng_t
  
  search_type = get_hdf5_type<u8>();
  is_equal = ( H5Tequal(datatype, search_type) > 0 );
  H5Tclose(search_type);
  if(is_equal) { return true; }
  
  search_type = get_hdf5_type<s8>();
  is_equal = ( H5Tequal(datatype, search_type) > 0 );
  H5Tclose(search_type);
  if(is_equal) { return true; }
  
  search_type = get_hdf5_type<u16>();
  is_equal = ( H5Tequal(datatype, search_type) > 0 );
  H5Tclose(search_type);
  if(is_equal) { return true; }
  
  search_type = get_hdf5_type<s16>();
  is_equal = ( H5Tequal(datatype, search_type) > 0 );
  H5Tclose(search_type);
  if(is_equal) { return true; }
  
  search_type = get_hdf5_type<u32>();
  is_equal = ( H5Tequal(datatype, search_type) > 0 );
  H5Tclose(search_type);
  if(is_equal) { return true; }
  
  search_type = get_hdf5_type<s32>();
  is_equal = ( H5Tequal(datatype, search_type) > 0 );
  H5Tclose(search_type);
  if(is_equal) { return true; }
  
  search_type = get_hdf5_type<u64>();
  is_equal = ( H5Tequal(datatype, search_type) > 0 );
  H5Tclose(search_type);
  if(is_equal) { return true; }
  
  search_type = get_hdf5_type<s64>();
  is_equal = ( H5Tequal(datatype, search_type) > 0 );
  H5Tclose(search_type);
  if(is_equal) { return true; }
  
  search_type = get_hdf5_type<ulng_t>();
  is_equal = ( H5Tequal(datatype, search_type) > 0 );
  H5Tclose(search_type);
  if(is_equal) { return true; }
  
  search_type = get_hdf5_type<slng_t>();
  is_equal = ( H5Tequal(datatype, search_type) > 0 );
  H5Tclose(search_type);
  if(is_equal) { return true; }
  
  return false;
  }



//! Auxiliary functions and structs for search_hdf5_file.
struct hdf5_search_info
  {
  const  std::vector<std::string>& names;
  int    num_dims;
  bool   exact;
  hid_t  best_match;
  size_t best_match_position; // Position of best match in names vector.
  };



inline
herr_t
hdf5_search_callback
  (
  hid_t             loc_id,
  const char*       name,
  const H5O_info_t* info,
  void*             operator_data  // hdf5_search_info
  )
  {
  hdf5_search_info* search_info = (hdf5_search_info*) operator_data;

  // We are looking for datasets.
  if(info->type == H5O_TYPE_DATASET)
    {
    // Check type of dataset to see if we could even load it.
    hid_t dataset  = H5Dopen(loc_id, name, H5P_DEFAULT);
    hid_t datatype = H5Dget_type(dataset);
    
    const bool is_supported = is_supported_arma_hdf5_type(datatype);
    
    H5Tclose(datatype);
    H5Dclose(dataset);
    
    if(is_supported == false)
      {
      // Forget about it and move on.
      return 0;
      }
    
    // Now we have to check against our set of names.
    // Only check names which could be better.
    for(size_t string_pos = 0; string_pos < search_info->best_match_position; ++string_pos)
      {
      // name is the full path (/path/to/dataset); names[string_pos] may be
      // "dataset", "/to/dataset", or "/path/to/dataset".
      // So if we count the number of forward slashes in names[string_pos],
      // and then simply take the last substring of name containing that number of slashes,
      // we can do the comparison.
      
      // Count the number of forward slashes in names[string_pos].
      uword name_count = 0;
      for(uword i = 0; i < search_info->names[string_pos].length(); ++i)
        {
        if((search_info->names[string_pos])[i] == '/') { ++name_count; }
        }

      // Count the number of forward slashes in the full name.
      uword count = 0;
      const std::string str = std::string(name);
      for(uword i = 0; i < str.length(); ++i)
        {
        if(str[i] == '/') { ++count; }
        }

      // Is the full string the same?
      if(str == search_info->names[string_pos])
        {
        // We found it exactly.
        hid_t match_candidate = H5Dopen(loc_id, name, H5P_DEFAULT);

        if(match_candidate < 0)
          {
          return -1;
          }

        // Ensure that the dataset is valid and of the correct dimensionality.
        hid_t filespace = H5Dget_space(match_candidate);
        int num_dims = H5Sget_simple_extent_ndims(filespace);
        
        if(num_dims <= search_info->num_dims)
          {
          // Valid dataset -- we'll keep it.
          // If we already have an existing match we have to close it.
          if(search_info->best_match != -1)
            {
            H5Dclose(search_info->best_match);
            }

          search_info->best_match_position = string_pos;
          search_info->best_match          = match_candidate;
          }
        
        H5Sclose(filespace);
        // There is no possibility of anything better, so terminate the search.
        return 1;
        }

      // If we are asking for more slashes than we have, this can't be a match.
      // Skip to below, where we decide whether or not to keep it anyway based
      // on the exactness condition of the search.
      if(count <= name_count)
        {
        size_t start_pos = (count == 0) ? 0 : std::string::npos;
        while(count > 0)
          {
          // Move pointer to previous slash.
          start_pos = str.rfind('/', start_pos);
          
          // Break if we've run out of slashes.
          if(start_pos == std::string::npos) { break; }
          
          --count;
          }

        // Now take the substring (this may end up being the full string).
        const std::string substring = str.substr(start_pos);

        // Are they the same?
        if(substring == search_info->names[string_pos])
          {
          // We have found the object; it must be better than our existing match.
          hid_t match_candidate = H5Dopen(loc_id, name, H5P_DEFAULT);


          // arma_check(match_candidate < 0, "Mat::load(): cannot open an HDF5 dataset");
          if(match_candidate < 0)
            {
            return -1;
            }
          
          
          // Ensure that the dataset is valid and of the correct dimensionality.
          hid_t filespace = H5Dget_space(match_candidate);
          int num_dims = H5Sget_simple_extent_ndims(filespace);
          
          if(num_dims <= search_info->num_dims)
            {
            // Valid dataset -- we'll keep it.
            // If we already have an existing match we have to close it.
            if(search_info->best_match != -1)
              {
              H5Dclose(search_info->best_match);
              }

            search_info->best_match_position = string_pos;
            search_info->best_match          = match_candidate;
            }
          
          H5Sclose(filespace);
          }
        }
      
      
      // If they are not the same, but we have not found anything and we don't need an exact match, take this.
      if((search_info->exact == false) && (search_info->best_match == -1))
        {
        hid_t match_candidate = H5Dopen(loc_id, name, H5P_DEFAULT);
        
        // arma_check(match_candidate < 0, "Mat::load(): cannot open an HDF5 dataset");
        if(match_candidate < 0)
          {
          return -1;
          }
        
        hid_t filespace = H5Dget_space(match_candidate);
        int num_dims = H5Sget_simple_extent_ndims(filespace);
        
        if(num_dims <= search_info->num_dims)
          {
          // Valid dataset -- we'll keep it.
          search_info->best_match = H5Dopen(loc_id, name, H5P_DEFAULT);
          }
        
        H5Sclose(filespace);
        }
      }
    }
  
  return 0;
  }



//! Search an HDF5 file for the given dataset names.  
//! If 'exact' is true, failure to find a dataset in the list of names means that -1 is returned.
//! If 'exact' is false and no datasets are found, -1 is returned.  
//! The number of dimensions is used to help prune down invalid datasets;
//! 2 dimensions is a matrix, 1 dimension is a vector, and 3 dimensions is a cube.  
//! If the number of dimensions in a dataset is less than or equal to num_dims, 
//! it will be considered -- for instance, a one-dimensional HDF5 vector can be loaded as a single-column matrix.
inline
hid_t
search_hdf5_file
  (
  const std::vector<std::string>& names,
  hid_t                           hdf5_file,
  int                             num_dims = 2,
  bool                            exact = false
  )
  {
  hdf5_search_info search_info = { names, num_dims, exact, -1, names.size() };
  
  // We'll use the H5Ovisit to track potential entries.
  herr_t status = H5Ovisit(hdf5_file, H5_INDEX_NAME, H5_ITER_NATIVE, hdf5_search_callback, void_ptr(&search_info));
  
  // Return the best match; it will be -1 if there was a problem.
  return (status < 0) ? -1 : search_info.best_match;
  }



//! Load an HDF5 matrix into an array of type specified by datatype,
//! then convert that into the desired array 'dest'.
//! This should only be called when eT is not the datatype.
template<typename eT>
inline
hid_t
load_and_convert_hdf5
  (
  eT   *dest,
  hid_t dataset,
  hid_t datatype,
  uword n_elem
  )
  {
  
  // We can't use nice template specializations here
  // as the determination of the type of 'datatype' must be done at runtime.
  // So we end up with this ugliness...
  hid_t search_type;
  
  bool is_equal;
  
  
  // u8
  search_type = get_hdf5_type<u8>();
  is_equal = (H5Tequal(datatype, search_type) > 0);
  H5Tclose(search_type);
  
  if(is_equal)
    {
    Col<u8> v(n_elem, arma_nozeros_indicator());
    hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
    arrayops::convert(dest, v.memptr(), n_elem);

    return status;
    }
  
  
  // s8
  search_type = get_hdf5_type<s8>();
  is_equal = (H5Tequal(datatype, search_type) > 0);
  H5Tclose(search_type);
  
  if(is_equal)
    {
    Col<s8> v(n_elem, arma_nozeros_indicator());
    hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
    arrayops::convert(dest, v.memptr(), n_elem);

    return status;
    }
  
  
  // u16
  search_type = get_hdf5_type<u16>();
  is_equal = (H5Tequal(datatype, search_type) > 0);
  H5Tclose(search_type);
  
  if(is_equal)
    {
    Col<u16> v(n_elem, arma_nozeros_indicator());
    hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
    arrayops::convert(dest, v.memptr(), n_elem);

    return status;
    }
  
  
  // s16
  search_type = get_hdf5_type<s16>();
  is_equal = (H5Tequal(datatype, search_type) > 0);
  H5Tclose(search_type);
  
  if(is_equal)
    {
    Col<s16> v(n_elem, arma_nozeros_indicator());
    hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
    arrayops::convert(dest, v.memptr(), n_elem);

    return status;
    }
  
  
  // u32
  search_type = get_hdf5_type<u32>();
  is_equal = (H5Tequal(datatype, search_type) > 0);
  H5Tclose(search_type);
  
  if(is_equal)
    {
    Col<u32> v(n_elem, arma_nozeros_indicator());
    hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
    arrayops::convert(dest, v.memptr(), n_elem);

    return status;
    }
  
  
  // s32
  search_type = get_hdf5_type<s32>();
  is_equal = (H5Tequal(datatype, search_type) > 0);
  H5Tclose(search_type);
  
  if(is_equal)
    {
    Col<s32> v(n_elem, arma_nozeros_indicator());
    hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
    arrayops::convert(dest, v.memptr(), n_elem);

    return status;
    }
  
  
  // u64
  search_type = get_hdf5_type<u64>();
  is_equal = (H5Tequal(datatype, search_type) > 0);
  H5Tclose(search_type);
  
  if(is_equal)
    {
    Col<u64> v(n_elem, arma_nozeros_indicator());
    hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
    arrayops::convert(dest, v.memptr(), n_elem);

    return status;
    }
  
  
  // s64
  search_type = get_hdf5_type<s64>();
  is_equal = (H5Tequal(datatype, search_type) > 0);
  H5Tclose(search_type);
  
  if(is_equal)
    {
    Col<s64> v(n_elem, arma_nozeros_indicator());
    hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
    arrayops::convert(dest, v.memptr(), n_elem);

    return status;
    }
  
  
  // ulng_t
  search_type = get_hdf5_type<ulng_t>();
  is_equal = (H5Tequal(datatype, search_type) > 0);
  H5Tclose(search_type);
  
  if(is_equal)
    {
    Col<ulng_t> v(n_elem, arma_nozeros_indicator());
    hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
    arrayops::convert(dest, v.memptr(), n_elem);

    return status;
    }
  
  
  // slng_t
  search_type = get_hdf5_type<slng_t>();
  is_equal = (H5Tequal(datatype, search_type) > 0);
  H5Tclose(search_type);
  
  if(is_equal)
    {
    Col<slng_t> v(n_elem, arma_nozeros_indicator());
    hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
    arrayops::convert(dest, v.memptr(), n_elem);

    return status;
    }
  
  
  // float
  search_type = get_hdf5_type<float>();
  is_equal = (H5Tequal(datatype, search_type) > 0);
  H5Tclose(search_type);
  
  if(is_equal)
    {
    Col<float> v(n_elem, arma_nozeros_indicator());
    hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
    arrayops::convert(dest, v.memptr(), n_elem);

    return status;
    }
  
  
  // double
  search_type = get_hdf5_type<double>();
  is_equal = (H5Tequal(datatype, search_type) > 0);
  H5Tclose(search_type);
  
  if(is_equal)
    {
    Col<double> v(n_elem, arma_nozeros_indicator());
    hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
    arrayops::convert(dest, v.memptr(), n_elem);

    return status;
    }
  
  
  // complex float
  search_type = get_hdf5_type< std::complex<float> >();
  is_equal = (H5Tequal(datatype, search_type) > 0);
  H5Tclose(search_type);
  
  if(is_equal)
    {
    if(is_cx<eT>::no)
      {
      return -1; // can't read complex data into non-complex matrix/cube
      }
    
    Col< std::complex<float> > v(n_elem, arma_nozeros_indicator());
    hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
    arrayops::convert_cx(dest, v.memptr(), n_elem);
    
    return status;
    }
  
  
  // complex double
  search_type = get_hdf5_type< std::complex<double> >();
  is_equal = (H5Tequal(datatype, search_type) > 0);
  H5Tclose(search_type);
  
  if(is_equal)
    {
    if(is_cx<eT>::no)
      {
      return -1; // can't read complex data into non-complex matrix/cube
      }
    
    Col< std::complex<double> > v(n_elem, arma_nozeros_indicator());
    hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
    arrayops::convert_cx(dest, v.memptr(), n_elem);
    
    return status;
    }
  
  
  return -1; // Failure.
  }



struct hdf5_suspend_printing_errors
  {
  #if (ARMA_WARN_LEVEL >= 3)
    
    inline
    hdf5_suspend_printing_errors() {}
    
  #else
    
    herr_t (*old_client_func)(hid_t, void*);
    void*    old_client_data;
    
    inline
    hdf5_suspend_printing_errors()
      {
      // Save old error handler.
      H5Eget_auto(H5E_DEFAULT, &old_client_func, &old_client_data);
      
      // Disable annoying HDF5 error messages.
      H5Eset_auto(H5E_DEFAULT, NULL, NULL);
      }
    
    inline
    ~hdf5_suspend_printing_errors()
      {
      H5Eset_auto(H5E_DEFAULT, old_client_func, old_client_data);
      }
    
  #endif
  };



}       // namespace hdf5_misc
#endif  // #if defined(ARMA_USE_HDF5)



//! @}