Chris@49
|
1 // Copyright (C) 2012 Ryan Curtin
|
Chris@49
|
2 // Copyright (C) 2012 Conrad Sanderson
|
Chris@49
|
3 // Copyright (C) 2013 Szabolcs Horvat
|
Chris@49
|
4 //
|
Chris@49
|
5 // This Source Code Form is subject to the terms of the Mozilla Public
|
Chris@49
|
6 // License, v. 2.0. If a copy of the MPL was not distributed with this
|
Chris@49
|
7 // file, You can obtain one at http://mozilla.org/MPL/2.0/.
|
Chris@49
|
8
|
Chris@49
|
9 //! \addtogroup hdf5_misc
|
Chris@49
|
10 //! @{
|
Chris@49
|
11
|
Chris@49
|
12
|
Chris@49
|
13 #if defined(ARMA_USE_HDF5)
|
Chris@49
|
14 namespace hdf5_misc
|
Chris@49
|
15 {
|
Chris@49
|
16
|
Chris@49
|
17
|
Chris@49
|
18 //! Given a certain type, find the corresponding HDF5 datatype. This can't be
|
Chris@49
|
19 //! done entirely at compile time, unfortunately, because the H5T_* macros
|
Chris@49
|
20 //! depend on function calls.
|
Chris@49
|
21 template< typename eT >
|
Chris@49
|
22 inline
|
Chris@49
|
23 hid_t
|
Chris@49
|
24 get_hdf5_type()
|
Chris@49
|
25 {
|
Chris@49
|
26 return -1; // Return invalid.
|
Chris@49
|
27 }
|
Chris@49
|
28
|
Chris@49
|
29
|
Chris@49
|
30
|
Chris@49
|
31 //! Specializations for each valid element type
|
Chris@49
|
32 //! (taken from all the possible typedefs of {u8, s8, ..., u64, s64} and the other native types.
|
Chris@49
|
33 //! We can't use the actual u8/s8 typedefs because their relations to the H5T_... types are unclear.
|
Chris@49
|
34 template<>
|
Chris@49
|
35 inline
|
Chris@49
|
36 hid_t
|
Chris@49
|
37 get_hdf5_type< unsigned char >()
|
Chris@49
|
38 {
|
Chris@49
|
39 return H5Tcopy(H5T_NATIVE_UCHAR);
|
Chris@49
|
40 }
|
Chris@49
|
41
|
Chris@49
|
42 template<>
|
Chris@49
|
43 inline
|
Chris@49
|
44 hid_t
|
Chris@49
|
45 get_hdf5_type< char >()
|
Chris@49
|
46 {
|
Chris@49
|
47 return H5Tcopy(H5T_NATIVE_CHAR);
|
Chris@49
|
48 }
|
Chris@49
|
49
|
Chris@49
|
50 template<>
|
Chris@49
|
51 inline
|
Chris@49
|
52 hid_t
|
Chris@49
|
53 get_hdf5_type< short >()
|
Chris@49
|
54 {
|
Chris@49
|
55 return H5Tcopy(H5T_NATIVE_SHORT);
|
Chris@49
|
56 }
|
Chris@49
|
57
|
Chris@49
|
58 template<>
|
Chris@49
|
59 inline
|
Chris@49
|
60 hid_t
|
Chris@49
|
61 get_hdf5_type< unsigned short >()
|
Chris@49
|
62 {
|
Chris@49
|
63 return H5Tcopy(H5T_NATIVE_USHORT);
|
Chris@49
|
64 }
|
Chris@49
|
65
|
Chris@49
|
66 template<>
|
Chris@49
|
67 inline
|
Chris@49
|
68 hid_t
|
Chris@49
|
69 get_hdf5_type< int >()
|
Chris@49
|
70 {
|
Chris@49
|
71 return H5Tcopy(H5T_NATIVE_INT);
|
Chris@49
|
72 }
|
Chris@49
|
73
|
Chris@49
|
74 template<>
|
Chris@49
|
75 inline
|
Chris@49
|
76 hid_t
|
Chris@49
|
77 get_hdf5_type< unsigned int >()
|
Chris@49
|
78 {
|
Chris@49
|
79 return H5Tcopy(H5T_NATIVE_UINT);
|
Chris@49
|
80 }
|
Chris@49
|
81
|
Chris@49
|
82 template<>
|
Chris@49
|
83 inline
|
Chris@49
|
84 hid_t
|
Chris@49
|
85 get_hdf5_type< long >()
|
Chris@49
|
86 {
|
Chris@49
|
87 return H5Tcopy(H5T_NATIVE_LONG);
|
Chris@49
|
88 }
|
Chris@49
|
89
|
Chris@49
|
90 template<>
|
Chris@49
|
91 inline
|
Chris@49
|
92 hid_t
|
Chris@49
|
93 get_hdf5_type< unsigned long >()
|
Chris@49
|
94 {
|
Chris@49
|
95 return H5Tcopy(H5T_NATIVE_ULONG);
|
Chris@49
|
96 }
|
Chris@49
|
97
|
Chris@49
|
98
|
Chris@49
|
99 #if defined(ARMA_USE_U64S64) && defined(ULLONG_MAX)
|
Chris@49
|
100 template<>
|
Chris@49
|
101 inline
|
Chris@49
|
102 hid_t
|
Chris@49
|
103 get_hdf5_type< long long >()
|
Chris@49
|
104 {
|
Chris@49
|
105 return H5Tcopy(H5T_NATIVE_LLONG);
|
Chris@49
|
106 }
|
Chris@49
|
107
|
Chris@49
|
108 template<>
|
Chris@49
|
109 inline
|
Chris@49
|
110 hid_t
|
Chris@49
|
111 get_hdf5_type< unsigned long long >()
|
Chris@49
|
112 {
|
Chris@49
|
113 return H5Tcopy(H5T_NATIVE_ULLONG);
|
Chris@49
|
114 }
|
Chris@49
|
115 #endif
|
Chris@49
|
116
|
Chris@49
|
117
|
Chris@49
|
118 template<>
|
Chris@49
|
119 inline
|
Chris@49
|
120 hid_t
|
Chris@49
|
121 get_hdf5_type< float >()
|
Chris@49
|
122 {
|
Chris@49
|
123 return H5Tcopy(H5T_NATIVE_FLOAT);
|
Chris@49
|
124 }
|
Chris@49
|
125
|
Chris@49
|
126 template<>
|
Chris@49
|
127 inline
|
Chris@49
|
128 hid_t
|
Chris@49
|
129 get_hdf5_type< double >()
|
Chris@49
|
130 {
|
Chris@49
|
131 return H5Tcopy(H5T_NATIVE_DOUBLE);
|
Chris@49
|
132 }
|
Chris@49
|
133
|
Chris@49
|
134
|
Chris@49
|
135
|
Chris@49
|
136 //! Utility hid_t since HOFFSET() won't work with std::complex.
|
Chris@49
|
137 template<typename eT>
|
Chris@49
|
138 struct hdf5_complex_t
|
Chris@49
|
139 {
|
Chris@49
|
140 eT real;
|
Chris@49
|
141 eT imag;
|
Chris@49
|
142 };
|
Chris@49
|
143
|
Chris@49
|
144
|
Chris@49
|
145
|
Chris@49
|
146 template<>
|
Chris@49
|
147 inline
|
Chris@49
|
148 hid_t
|
Chris@49
|
149 get_hdf5_type< std::complex<float> >()
|
Chris@49
|
150 {
|
Chris@49
|
151 hid_t type = H5Tcreate(H5T_COMPOUND, sizeof(hdf5_complex_t<float>));
|
Chris@49
|
152
|
Chris@49
|
153 H5Tinsert(type, "real", HOFFSET(hdf5_complex_t<float>, real), H5T_NATIVE_FLOAT);
|
Chris@49
|
154 H5Tinsert(type, "imag", HOFFSET(hdf5_complex_t<float>, imag), H5T_NATIVE_FLOAT);
|
Chris@49
|
155
|
Chris@49
|
156 return type;
|
Chris@49
|
157 }
|
Chris@49
|
158
|
Chris@49
|
159
|
Chris@49
|
160
|
Chris@49
|
161 template<>
|
Chris@49
|
162 inline
|
Chris@49
|
163 hid_t
|
Chris@49
|
164 get_hdf5_type< std::complex<double> >()
|
Chris@49
|
165 {
|
Chris@49
|
166 hid_t type = H5Tcreate(H5T_COMPOUND, sizeof(hdf5_complex_t<double>));
|
Chris@49
|
167
|
Chris@49
|
168 H5Tinsert(type, "real", HOFFSET(hdf5_complex_t<double>, real), H5T_NATIVE_DOUBLE);
|
Chris@49
|
169 H5Tinsert(type, "imag", HOFFSET(hdf5_complex_t<double>, imag), H5T_NATIVE_DOUBLE);
|
Chris@49
|
170
|
Chris@49
|
171 return type;
|
Chris@49
|
172 }
|
Chris@49
|
173
|
Chris@49
|
174
|
Chris@49
|
175
|
Chris@49
|
176 // Compare datatype against all supported types.
|
Chris@49
|
177 inline
|
Chris@49
|
178 bool
|
Chris@49
|
179 is_supported_arma_hdf5_type(hid_t datatype)
|
Chris@49
|
180 {
|
Chris@49
|
181 hid_t search_type;
|
Chris@49
|
182
|
Chris@49
|
183 bool is_equal;
|
Chris@49
|
184
|
Chris@49
|
185
|
Chris@49
|
186 // start with most likely used types: double, complex<double>, float, complex<float>
|
Chris@49
|
187
|
Chris@49
|
188 search_type = get_hdf5_type<double>();
|
Chris@49
|
189 is_equal = ( H5Tequal(datatype, search_type) > 0 );
|
Chris@49
|
190 H5Tclose(search_type);
|
Chris@49
|
191 if (is_equal) { return true; }
|
Chris@49
|
192
|
Chris@49
|
193 search_type = get_hdf5_type< std::complex<double> >();
|
Chris@49
|
194 is_equal = ( H5Tequal(datatype, search_type) > 0 );
|
Chris@49
|
195 H5Tclose(search_type);
|
Chris@49
|
196 if (is_equal) { return true; }
|
Chris@49
|
197
|
Chris@49
|
198 search_type = get_hdf5_type<float>();
|
Chris@49
|
199 is_equal = ( H5Tequal(datatype, search_type) > 0 );
|
Chris@49
|
200 H5Tclose(search_type);
|
Chris@49
|
201 if (is_equal) { return true; }
|
Chris@49
|
202
|
Chris@49
|
203 search_type = get_hdf5_type< std::complex<float> >();
|
Chris@49
|
204 is_equal = ( H5Tequal(datatype, search_type) > 0 );
|
Chris@49
|
205 H5Tclose(search_type);
|
Chris@49
|
206 if (is_equal) { return true; }
|
Chris@49
|
207
|
Chris@49
|
208
|
Chris@49
|
209 // remaining supported types: u8, s8, u16, s16, u32, s32, u64, s64, ulng_t, slng_t
|
Chris@49
|
210
|
Chris@49
|
211 search_type = get_hdf5_type<u8>();
|
Chris@49
|
212 is_equal = ( H5Tequal(datatype, search_type) > 0 );
|
Chris@49
|
213 H5Tclose(search_type);
|
Chris@49
|
214 if (is_equal) { return true; }
|
Chris@49
|
215
|
Chris@49
|
216 search_type = get_hdf5_type<s8>();
|
Chris@49
|
217 is_equal = ( H5Tequal(datatype, search_type) > 0 );
|
Chris@49
|
218 H5Tclose(search_type);
|
Chris@49
|
219 if (is_equal) { return true; }
|
Chris@49
|
220
|
Chris@49
|
221 search_type = get_hdf5_type<u16>();
|
Chris@49
|
222 is_equal = ( H5Tequal(datatype, search_type) > 0 );
|
Chris@49
|
223 H5Tclose(search_type);
|
Chris@49
|
224 if (is_equal) { return true; }
|
Chris@49
|
225
|
Chris@49
|
226 search_type = get_hdf5_type<s16>();
|
Chris@49
|
227 is_equal = ( H5Tequal(datatype, search_type) > 0 );
|
Chris@49
|
228 H5Tclose(search_type);
|
Chris@49
|
229 if (is_equal) { return true; }
|
Chris@49
|
230
|
Chris@49
|
231 search_type = get_hdf5_type<u32>();
|
Chris@49
|
232 is_equal = ( H5Tequal(datatype, search_type) > 0 );
|
Chris@49
|
233 H5Tclose(search_type);
|
Chris@49
|
234 if (is_equal) { return true; }
|
Chris@49
|
235
|
Chris@49
|
236 search_type = get_hdf5_type<s32>();
|
Chris@49
|
237 is_equal = ( H5Tequal(datatype, search_type) > 0 );
|
Chris@49
|
238 H5Tclose(search_type);
|
Chris@49
|
239 if (is_equal) { return true; }
|
Chris@49
|
240
|
Chris@49
|
241 #if defined(ARMA_USE_U64S64)
|
Chris@49
|
242 {
|
Chris@49
|
243 search_type = get_hdf5_type<u64>();
|
Chris@49
|
244 is_equal = ( H5Tequal(datatype, search_type) > 0 );
|
Chris@49
|
245 H5Tclose(search_type);
|
Chris@49
|
246 if (is_equal) { return true; }
|
Chris@49
|
247
|
Chris@49
|
248 search_type = get_hdf5_type<s64>();
|
Chris@49
|
249 is_equal = ( H5Tequal(datatype, search_type) > 0 );
|
Chris@49
|
250 H5Tclose(search_type);
|
Chris@49
|
251 if (is_equal) { return true; }
|
Chris@49
|
252 }
|
Chris@49
|
253 #endif
|
Chris@49
|
254
|
Chris@49
|
255 #if defined(ARMA_ALLOW_LONG)
|
Chris@49
|
256 {
|
Chris@49
|
257 search_type = get_hdf5_type<ulng_t>();
|
Chris@49
|
258 is_equal = ( H5Tequal(datatype, search_type) > 0 );
|
Chris@49
|
259 H5Tclose(search_type);
|
Chris@49
|
260 if (is_equal) { return true; }
|
Chris@49
|
261
|
Chris@49
|
262 search_type = get_hdf5_type<slng_t>();
|
Chris@49
|
263 is_equal = ( H5Tequal(datatype, search_type) > 0 );
|
Chris@49
|
264 H5Tclose(search_type);
|
Chris@49
|
265 if (is_equal) { return true; }
|
Chris@49
|
266 }
|
Chris@49
|
267 #endif
|
Chris@49
|
268
|
Chris@49
|
269 return false;
|
Chris@49
|
270 }
|
Chris@49
|
271
|
Chris@49
|
272
|
Chris@49
|
273
|
Chris@49
|
274 //! Auxiliary functions and structs for search_hdf5_file.
|
Chris@49
|
275 struct hdf5_search_info
|
Chris@49
|
276 {
|
Chris@49
|
277 const std::vector<std::string>& names;
|
Chris@49
|
278 int num_dims;
|
Chris@49
|
279 bool exact;
|
Chris@49
|
280 hid_t best_match;
|
Chris@49
|
281 size_t best_match_position; // Position of best match in names vector.
|
Chris@49
|
282 };
|
Chris@49
|
283
|
Chris@49
|
284
|
Chris@49
|
285
|
Chris@49
|
286 inline
|
Chris@49
|
287 herr_t
|
Chris@49
|
288 hdf5_search_callback
|
Chris@49
|
289 (
|
Chris@49
|
290 hid_t loc_id,
|
Chris@49
|
291 const char* name,
|
Chris@49
|
292 const H5O_info_t* info,
|
Chris@49
|
293 void* operator_data // hdf5_search_info
|
Chris@49
|
294 )
|
Chris@49
|
295 {
|
Chris@49
|
296 hdf5_search_info* search_info = (hdf5_search_info*) operator_data;
|
Chris@49
|
297
|
Chris@49
|
298 // We are looking for datasets.
|
Chris@49
|
299 if (info->type == H5O_TYPE_DATASET)
|
Chris@49
|
300 {
|
Chris@49
|
301 // Check type of dataset to see if we could even load it.
|
Chris@49
|
302 hid_t dataset = H5Dopen(loc_id, name, H5P_DEFAULT);
|
Chris@49
|
303 hid_t datatype = H5Dget_type(dataset);
|
Chris@49
|
304
|
Chris@49
|
305 const bool is_supported = is_supported_arma_hdf5_type(datatype);
|
Chris@49
|
306
|
Chris@49
|
307 H5Tclose(datatype);
|
Chris@49
|
308 H5Dclose(dataset);
|
Chris@49
|
309
|
Chris@49
|
310 if(is_supported == false)
|
Chris@49
|
311 {
|
Chris@49
|
312 // Forget about it and move on.
|
Chris@49
|
313 return 0;
|
Chris@49
|
314 }
|
Chris@49
|
315
|
Chris@49
|
316 // Now we have to check against our set of names.
|
Chris@49
|
317 // Only check names which could be better.
|
Chris@49
|
318 for (size_t string_pos = 0; string_pos < search_info->best_match_position; ++string_pos)
|
Chris@49
|
319 {
|
Chris@49
|
320 // name is the full path (/path/to/dataset); names[string_pos] may be
|
Chris@49
|
321 // "dataset", "/to/dataset", or "/path/to/dataset".
|
Chris@49
|
322 // So if we count the number of forward slashes in names[string_pos],
|
Chris@49
|
323 // and then simply take the last substring of name containing that number of slashes,
|
Chris@49
|
324 // we can do the comparison.
|
Chris@49
|
325
|
Chris@49
|
326 // Count the number of forward slashes in names[string_pos].
|
Chris@49
|
327 uword count = 0;
|
Chris@49
|
328 for (uword i = 0; i < search_info->names[string_pos].length(); ++i)
|
Chris@49
|
329 {
|
Chris@49
|
330 if ((search_info->names[string_pos])[i] == '/') { ++count; }
|
Chris@49
|
331 }
|
Chris@49
|
332
|
Chris@49
|
333 // Count the number of forward slashes in the full name.
|
Chris@49
|
334 uword name_count = 0;
|
Chris@49
|
335 const std::string str = std::string(name);
|
Chris@49
|
336 for (uword i = 0; i < str.length(); ++i)
|
Chris@49
|
337 {
|
Chris@49
|
338 if (str[i] == '/') { ++count; }
|
Chris@49
|
339 }
|
Chris@49
|
340
|
Chris@49
|
341 // If we are asking for more slashes than we have, this can't be a match.
|
Chris@49
|
342 // Skip to below, where we decide whether or not to keep it anyway based
|
Chris@49
|
343 // on the exactness condition of the search.
|
Chris@49
|
344 if (count <= name_count)
|
Chris@49
|
345 {
|
Chris@49
|
346 size_t start_pos = (count == 0) ? 0 : std::string::npos;
|
Chris@49
|
347 while (count > 0)
|
Chris@49
|
348 {
|
Chris@49
|
349 // Move pointer to previous slash.
|
Chris@49
|
350 start_pos = str.rfind('/', start_pos);
|
Chris@49
|
351
|
Chris@49
|
352 // Break if we've run out of slashes.
|
Chris@49
|
353 if (start_pos == std::string::npos) { break; }
|
Chris@49
|
354
|
Chris@49
|
355 --count;
|
Chris@49
|
356 }
|
Chris@49
|
357
|
Chris@49
|
358 // Now take the substring (this may end up being the full string).
|
Chris@49
|
359 const std::string substring = str.substr(start_pos);
|
Chris@49
|
360
|
Chris@49
|
361 // Are they the same?
|
Chris@49
|
362 if (substring == search_info->names[string_pos])
|
Chris@49
|
363 {
|
Chris@49
|
364 // We have found the object; it must be better than our existing match.
|
Chris@49
|
365 hid_t match_candidate = H5Dopen(loc_id, name, H5P_DEFAULT);
|
Chris@49
|
366
|
Chris@49
|
367
|
Chris@49
|
368 // arma_check(match_candidate < 0, "Mat::load(): cannot open an HDF5 dataset");
|
Chris@49
|
369 if(match_candidate < 0)
|
Chris@49
|
370 {
|
Chris@49
|
371 return -1;
|
Chris@49
|
372 }
|
Chris@49
|
373
|
Chris@49
|
374
|
Chris@49
|
375 // Ensure that the dataset is valid and of the correct dimensionality.
|
Chris@49
|
376 hid_t filespace = H5Dget_space(match_candidate);
|
Chris@49
|
377 int num_dims = H5Sget_simple_extent_ndims(filespace);
|
Chris@49
|
378
|
Chris@49
|
379 if (num_dims <= search_info->num_dims)
|
Chris@49
|
380 {
|
Chris@49
|
381 // Valid dataset -- we'll keep it.
|
Chris@49
|
382 // If we already have an existing match we have to close it.
|
Chris@49
|
383 if (search_info->best_match != -1)
|
Chris@49
|
384 {
|
Chris@49
|
385 H5Dclose(search_info->best_match);
|
Chris@49
|
386 }
|
Chris@49
|
387
|
Chris@49
|
388 search_info->best_match_position = string_pos;
|
Chris@49
|
389 search_info->best_match = match_candidate;
|
Chris@49
|
390 }
|
Chris@49
|
391
|
Chris@49
|
392 H5Sclose(filespace);
|
Chris@49
|
393 }
|
Chris@49
|
394 }
|
Chris@49
|
395
|
Chris@49
|
396
|
Chris@49
|
397 // If they are not the same, but we have not found anything and we don't need an exact match, take this.
|
Chris@49
|
398 if ((search_info->exact == false) && (search_info->best_match == -1))
|
Chris@49
|
399 {
|
Chris@49
|
400 hid_t match_candidate = H5Dopen(loc_id, name, H5P_DEFAULT);
|
Chris@49
|
401
|
Chris@49
|
402 // arma_check(match_candidate < 0, "Mat::load(): cannot open an HDF5 dataset");
|
Chris@49
|
403 if(match_candidate < 0)
|
Chris@49
|
404 {
|
Chris@49
|
405 return -1;
|
Chris@49
|
406 }
|
Chris@49
|
407
|
Chris@49
|
408 hid_t filespace = H5Dget_space(match_candidate);
|
Chris@49
|
409 int num_dims = H5Sget_simple_extent_ndims(filespace);
|
Chris@49
|
410
|
Chris@49
|
411 if (num_dims <= search_info->num_dims)
|
Chris@49
|
412 {
|
Chris@49
|
413 // Valid dataset -- we'll keep it.
|
Chris@49
|
414 search_info->best_match = H5Dopen(loc_id, name, H5P_DEFAULT);
|
Chris@49
|
415 }
|
Chris@49
|
416
|
Chris@49
|
417 H5Sclose(filespace);
|
Chris@49
|
418 }
|
Chris@49
|
419 }
|
Chris@49
|
420 }
|
Chris@49
|
421
|
Chris@49
|
422 return 0;
|
Chris@49
|
423 }
|
Chris@49
|
424
|
Chris@49
|
425
|
Chris@49
|
426
|
Chris@49
|
427 //! Search an HDF5 file for the given dataset names.
|
Chris@49
|
428 //! If 'exact' is true, failure to find a dataset in the list of names means that -1 is returned.
|
Chris@49
|
429 //! If 'exact' is false and no datasets are found, -1 is returned.
|
Chris@49
|
430 //! The number of dimensions is used to help prune down invalid datasets;
|
Chris@49
|
431 //! 2 dimensions is a matrix, 1 dimension is a vector, and 3 dimensions is a cube.
|
Chris@49
|
432 //! If the number of dimensions in a dataset is less than or equal to num_dims,
|
Chris@49
|
433 //! it will be considered -- for instance, a one-dimensional HDF5 vector can be loaded as a single-column matrix.
|
Chris@49
|
434 inline
|
Chris@49
|
435 hid_t
|
Chris@49
|
436 search_hdf5_file
|
Chris@49
|
437 (
|
Chris@49
|
438 const std::vector<std::string>& names,
|
Chris@49
|
439 hid_t hdf5_file,
|
Chris@49
|
440 int num_dims = 2,
|
Chris@49
|
441 bool exact = false
|
Chris@49
|
442 )
|
Chris@49
|
443 {
|
Chris@49
|
444 hdf5_search_info search_info = { names, num_dims, exact, -1, names.size() };
|
Chris@49
|
445
|
Chris@49
|
446 // We'll use the H5Ovisit to track potential entries.
|
Chris@49
|
447 herr_t status = H5Ovisit(hdf5_file, H5_INDEX_NAME, H5_ITER_NATIVE, hdf5_search_callback, void_ptr(&search_info));
|
Chris@49
|
448
|
Chris@49
|
449 // Return the best match; it will be -1 if there was a problem.
|
Chris@49
|
450 return (status < 0) ? -1 : search_info.best_match;
|
Chris@49
|
451 }
|
Chris@49
|
452
|
Chris@49
|
453
|
Chris@49
|
454
|
Chris@49
|
455 //! Load an HDF5 matrix into an array of type specified by datatype,
|
Chris@49
|
456 //! then convert that into the desired array 'dest'.
|
Chris@49
|
457 //! This should only be called when eT is not the datatype.
|
Chris@49
|
458 template<typename eT>
|
Chris@49
|
459 inline
|
Chris@49
|
460 hid_t
|
Chris@49
|
461 load_and_convert_hdf5
|
Chris@49
|
462 (
|
Chris@49
|
463 eT *dest,
|
Chris@49
|
464 hid_t dataset,
|
Chris@49
|
465 hid_t datatype,
|
Chris@49
|
466 uword n_elem
|
Chris@49
|
467 )
|
Chris@49
|
468 {
|
Chris@49
|
469
|
Chris@49
|
470 // We can't use nice template specializations here
|
Chris@49
|
471 // as the determination of the type of 'datatype' must be done at runtime.
|
Chris@49
|
472 // So we end up with this ugliness...
|
Chris@49
|
473 hid_t search_type;
|
Chris@49
|
474
|
Chris@49
|
475 bool is_equal;
|
Chris@49
|
476
|
Chris@49
|
477
|
Chris@49
|
478 // u8
|
Chris@49
|
479 search_type = get_hdf5_type<u8>();
|
Chris@49
|
480 is_equal = (H5Tequal(datatype, search_type) > 0);
|
Chris@49
|
481 H5Tclose(search_type);
|
Chris@49
|
482
|
Chris@49
|
483 if(is_equal)
|
Chris@49
|
484 {
|
Chris@49
|
485 Col<u8> v(n_elem);
|
Chris@49
|
486 hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
|
Chris@49
|
487 arrayops::convert(dest, v.memptr(), n_elem);
|
Chris@49
|
488
|
Chris@49
|
489 return status;
|
Chris@49
|
490 }
|
Chris@49
|
491
|
Chris@49
|
492
|
Chris@49
|
493 // s8
|
Chris@49
|
494 search_type = get_hdf5_type<s8>();
|
Chris@49
|
495 is_equal = (H5Tequal(datatype, search_type) > 0);
|
Chris@49
|
496 H5Tclose(search_type);
|
Chris@49
|
497
|
Chris@49
|
498 if(is_equal)
|
Chris@49
|
499 {
|
Chris@49
|
500 Col<s8> v(n_elem);
|
Chris@49
|
501 hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
|
Chris@49
|
502 arrayops::convert(dest, v.memptr(), n_elem);
|
Chris@49
|
503
|
Chris@49
|
504 return status;
|
Chris@49
|
505 }
|
Chris@49
|
506
|
Chris@49
|
507
|
Chris@49
|
508 // u16
|
Chris@49
|
509 search_type = get_hdf5_type<u16>();
|
Chris@49
|
510 is_equal = (H5Tequal(datatype, search_type) > 0);
|
Chris@49
|
511 H5Tclose(search_type);
|
Chris@49
|
512
|
Chris@49
|
513 if(is_equal)
|
Chris@49
|
514 {
|
Chris@49
|
515 Col<u16> v(n_elem);
|
Chris@49
|
516 hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
|
Chris@49
|
517 arrayops::convert(dest, v.memptr(), n_elem);
|
Chris@49
|
518
|
Chris@49
|
519 return status;
|
Chris@49
|
520 }
|
Chris@49
|
521
|
Chris@49
|
522
|
Chris@49
|
523 // s16
|
Chris@49
|
524 search_type = get_hdf5_type<s16>();
|
Chris@49
|
525 is_equal = (H5Tequal(datatype, search_type) > 0);
|
Chris@49
|
526 H5Tclose(search_type);
|
Chris@49
|
527
|
Chris@49
|
528 if(is_equal)
|
Chris@49
|
529 {
|
Chris@49
|
530 Col<s16> v(n_elem);
|
Chris@49
|
531 hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
|
Chris@49
|
532 arrayops::convert(dest, v.memptr(), n_elem);
|
Chris@49
|
533
|
Chris@49
|
534 return status;
|
Chris@49
|
535 }
|
Chris@49
|
536
|
Chris@49
|
537
|
Chris@49
|
538 // u32
|
Chris@49
|
539 search_type = get_hdf5_type<u32>();
|
Chris@49
|
540 is_equal = (H5Tequal(datatype, search_type) > 0);
|
Chris@49
|
541 H5Tclose(search_type);
|
Chris@49
|
542
|
Chris@49
|
543 if(is_equal)
|
Chris@49
|
544 {
|
Chris@49
|
545 Col<u32> v(n_elem);
|
Chris@49
|
546 hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
|
Chris@49
|
547 arrayops::convert(dest, v.memptr(), n_elem);
|
Chris@49
|
548
|
Chris@49
|
549 return status;
|
Chris@49
|
550 }
|
Chris@49
|
551
|
Chris@49
|
552
|
Chris@49
|
553 // s32
|
Chris@49
|
554 search_type = get_hdf5_type<s32>();
|
Chris@49
|
555 is_equal = (H5Tequal(datatype, search_type) > 0);
|
Chris@49
|
556 H5Tclose(search_type);
|
Chris@49
|
557
|
Chris@49
|
558 if(is_equal)
|
Chris@49
|
559 {
|
Chris@49
|
560 Col<s32> v(n_elem);
|
Chris@49
|
561 hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
|
Chris@49
|
562 arrayops::convert(dest, v.memptr(), n_elem);
|
Chris@49
|
563
|
Chris@49
|
564 return status;
|
Chris@49
|
565 }
|
Chris@49
|
566
|
Chris@49
|
567
|
Chris@49
|
568 #if defined(ARMA_USE_U64S64)
|
Chris@49
|
569 {
|
Chris@49
|
570 // u64
|
Chris@49
|
571 search_type = get_hdf5_type<u64>();
|
Chris@49
|
572 is_equal = (H5Tequal(datatype, search_type) > 0);
|
Chris@49
|
573 H5Tclose(search_type);
|
Chris@49
|
574
|
Chris@49
|
575 if(is_equal)
|
Chris@49
|
576 {
|
Chris@49
|
577 Col<u64> v(n_elem);
|
Chris@49
|
578 hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
|
Chris@49
|
579 arrayops::convert(dest, v.memptr(), n_elem);
|
Chris@49
|
580
|
Chris@49
|
581 return status;
|
Chris@49
|
582 }
|
Chris@49
|
583
|
Chris@49
|
584
|
Chris@49
|
585 // s64
|
Chris@49
|
586 search_type = get_hdf5_type<s64>();
|
Chris@49
|
587 is_equal = (H5Tequal(datatype, search_type) > 0);
|
Chris@49
|
588 H5Tclose(search_type);
|
Chris@49
|
589
|
Chris@49
|
590 if(is_equal)
|
Chris@49
|
591 {
|
Chris@49
|
592 Col<s64> v(n_elem);
|
Chris@49
|
593 hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
|
Chris@49
|
594 arrayops::convert(dest, v.memptr(), n_elem);
|
Chris@49
|
595
|
Chris@49
|
596 return status;
|
Chris@49
|
597 }
|
Chris@49
|
598 }
|
Chris@49
|
599 #endif
|
Chris@49
|
600
|
Chris@49
|
601
|
Chris@49
|
602 #if defined(ARMA_ALLOW_LONG)
|
Chris@49
|
603 {
|
Chris@49
|
604 // ulng_t
|
Chris@49
|
605 search_type = get_hdf5_type<ulng_t>();
|
Chris@49
|
606 is_equal = (H5Tequal(datatype, search_type) > 0);
|
Chris@49
|
607 H5Tclose(search_type);
|
Chris@49
|
608
|
Chris@49
|
609 if(is_equal)
|
Chris@49
|
610 {
|
Chris@49
|
611 Col<ulng_t> v(n_elem);
|
Chris@49
|
612 hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
|
Chris@49
|
613 arrayops::convert(dest, v.memptr(), n_elem);
|
Chris@49
|
614
|
Chris@49
|
615 return status;
|
Chris@49
|
616 }
|
Chris@49
|
617
|
Chris@49
|
618
|
Chris@49
|
619 // slng_t
|
Chris@49
|
620 search_type = get_hdf5_type<slng_t>();
|
Chris@49
|
621 is_equal = (H5Tequal(datatype, search_type) > 0);
|
Chris@49
|
622 H5Tclose(search_type);
|
Chris@49
|
623
|
Chris@49
|
624 if(is_equal)
|
Chris@49
|
625 {
|
Chris@49
|
626 Col<slng_t> v(n_elem);
|
Chris@49
|
627 hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
|
Chris@49
|
628 arrayops::convert(dest, v.memptr(), n_elem);
|
Chris@49
|
629
|
Chris@49
|
630 return status;
|
Chris@49
|
631 }
|
Chris@49
|
632 }
|
Chris@49
|
633 #endif
|
Chris@49
|
634
|
Chris@49
|
635
|
Chris@49
|
636 // float
|
Chris@49
|
637 search_type = get_hdf5_type<float>();
|
Chris@49
|
638 is_equal = (H5Tequal(datatype, search_type) > 0);
|
Chris@49
|
639 H5Tclose(search_type);
|
Chris@49
|
640
|
Chris@49
|
641 if(is_equal)
|
Chris@49
|
642 {
|
Chris@49
|
643 Col<float> v(n_elem);
|
Chris@49
|
644 hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
|
Chris@49
|
645 arrayops::convert(dest, v.memptr(), n_elem);
|
Chris@49
|
646
|
Chris@49
|
647 return status;
|
Chris@49
|
648 }
|
Chris@49
|
649
|
Chris@49
|
650
|
Chris@49
|
651 // double
|
Chris@49
|
652 search_type = get_hdf5_type<double>();
|
Chris@49
|
653 is_equal = (H5Tequal(datatype, search_type) > 0);
|
Chris@49
|
654 H5Tclose(search_type);
|
Chris@49
|
655
|
Chris@49
|
656 if(is_equal)
|
Chris@49
|
657 {
|
Chris@49
|
658 Col<double> v(n_elem);
|
Chris@49
|
659 hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
|
Chris@49
|
660 arrayops::convert(dest, v.memptr(), n_elem);
|
Chris@49
|
661
|
Chris@49
|
662 return status;
|
Chris@49
|
663 }
|
Chris@49
|
664
|
Chris@49
|
665
|
Chris@49
|
666 // complex float
|
Chris@49
|
667 search_type = get_hdf5_type< std::complex<float> >();
|
Chris@49
|
668 is_equal = (H5Tequal(datatype, search_type) > 0);
|
Chris@49
|
669 H5Tclose(search_type);
|
Chris@49
|
670
|
Chris@49
|
671 if(is_equal)
|
Chris@49
|
672 {
|
Chris@49
|
673 if(is_complex<eT>::value == false)
|
Chris@49
|
674 {
|
Chris@49
|
675 return -1; // can't read complex data into non-complex matrix/cube
|
Chris@49
|
676 }
|
Chris@49
|
677
|
Chris@49
|
678 Col< std::complex<float> > v(n_elem);
|
Chris@49
|
679 hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
|
Chris@49
|
680 arrayops::convert_cx(dest, v.memptr(), n_elem);
|
Chris@49
|
681
|
Chris@49
|
682 return status;
|
Chris@49
|
683 }
|
Chris@49
|
684
|
Chris@49
|
685
|
Chris@49
|
686 // complex double
|
Chris@49
|
687 search_type = get_hdf5_type< std::complex<double> >();
|
Chris@49
|
688 is_equal = (H5Tequal(datatype, search_type) > 0);
|
Chris@49
|
689 H5Tclose(search_type);
|
Chris@49
|
690
|
Chris@49
|
691 if(is_equal)
|
Chris@49
|
692 {
|
Chris@49
|
693 if(is_complex<eT>::value == false)
|
Chris@49
|
694 {
|
Chris@49
|
695 return -1; // can't read complex data into non-complex matrix/cube
|
Chris@49
|
696 }
|
Chris@49
|
697
|
Chris@49
|
698 Col< std::complex<double> > v(n_elem);
|
Chris@49
|
699 hid_t status = H5Dread(dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, void_ptr(v.memptr()));
|
Chris@49
|
700 arrayops::convert_cx(dest, v.memptr(), n_elem);
|
Chris@49
|
701
|
Chris@49
|
702 return status;
|
Chris@49
|
703 }
|
Chris@49
|
704
|
Chris@49
|
705
|
Chris@49
|
706 return -1; // Failure.
|
Chris@49
|
707 }
|
Chris@49
|
708
|
Chris@49
|
709
|
Chris@49
|
710
|
Chris@49
|
711 } // namespace hdf5_misc
|
Chris@49
|
712 #endif // #if defined(ARMA_USE_HDF5)
|
Chris@49
|
713
|
Chris@49
|
714
|
Chris@49
|
715
|
Chris@49
|
716 //! @}
|