forked from LLNL/libROM
-
Notifications
You must be signed in to change notification settings - Fork 0
/
HDFDatabase.C
393 lines (335 loc) · 9.26 KB
/
HDFDatabase.C
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
/******************************************************************************
*
* Copyright (c) 2013-2019, Lawrence Livermore National Security, LLC
* and other libROM project developers. See the top-level COPYRIGHT
* file for details.
*
* SPDX-License-Identifier: (Apache-2.0 OR MIT)
*
*****************************************************************************/
// Description: The concrete database implementation using HDF5.
#include "HDFDatabase.h"
#include "Utilities.h"
namespace CAROM {
const int HDFDatabase::KEY_DOUBLE_ARRAY = 0;
const int HDFDatabase::KEY_INT_ARRAY = 1;
HDFDatabase::HDFDatabase() :
d_is_file(false),
d_file_id(-1),
d_group_id(-1)
{
}
HDFDatabase::~HDFDatabase()
{
if (d_is_file) {
close();
}
if (d_group_id != -1) {
herr_t errf = H5Gclose(d_group_id);
CAROM_ASSERT(errf >= 0);
#ifndef DEBUG_CHECK_ASSERTIONS
CAROM_NULL_USE(errf);
#endif
}
}
bool
HDFDatabase::create(
const std::string& file_name)
{
CAROM_ASSERT(!file_name.empty());
hid_t file_id = H5Fcreate(file_name.c_str(),
H5F_ACC_TRUNC,
H5P_DEFAULT,
H5P_DEFAULT);
bool result = file_id >= 0;
CAROM_ASSERT(result);
d_is_file = true;
d_file_id = file_id;
d_group_id = file_id;
return result;
}
bool
HDFDatabase::open(
const std::string& file_name)
{
CAROM_ASSERT(!file_name.empty());
hid_t file_id = H5Fopen(file_name.c_str(),
H5F_ACC_RDONLY,
H5P_DEFAULT);
bool result = file_id >= 0;
CAROM_ASSERT(result);
d_is_file = true;
d_file_id = file_id;
d_group_id = file_id;
return result;
}
bool
HDFDatabase::close()
{
herr_t errf = 0;
if (d_is_file) {
errf = H5Fclose(d_file_id);
CAROM_ASSERT(errf >= 0);
if (d_group_id == d_file_id) {
d_group_id = -1;
}
d_file_id = -1;
d_is_file = false;
}
return errf >= 0;
}
void
HDFDatabase::putIntegerArray(
const std::string& key,
const int* const data,
int nelements)
{
CAROM_ASSERT(!key.empty());
CAROM_ASSERT(data != 0);
CAROM_ASSERT(nelements > 0);
hsize_t dim[] = { static_cast<hsize_t>(nelements) };
hid_t space = H5Screate_simple(1, dim, 0);
CAROM_ASSERT(space >= 0);
#if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
hid_t dataset = H5Dcreate(d_group_id,
key.c_str(),
H5T_STD_I32BE,
space,
H5P_DEFAULT,
H5P_DEFAULT,
H5P_DEFAULT);
#else
hid_t dataset = H5Dcreate(d_group_id,
key.c_str(),
H5T_STD_I32BE,
space,
H5P_DEFAULT);
#endif
CAROM_ASSERT(dataset >= 0);
herr_t errf = H5Dwrite(dataset,
H5T_NATIVE_INT,
H5S_ALL,
H5S_ALL,
H5P_DEFAULT,
data);
CAROM_ASSERT(errf >= 0);
// Write attribute so we know what kind of data this is.
writeAttribute(KEY_INT_ARRAY, dataset);
errf = H5Sclose(space);
CAROM_ASSERT(errf >= 0);
errf = H5Dclose(dataset);
CAROM_ASSERT(errf >= 0);
#ifndef DEBUG_CHECK_ASSERTIONS
CAROM_NULL_USE(errf);
#endif
}
void
HDFDatabase::putDoubleArray(
const std::string& key,
const double* const data,
int nelements)
{
CAROM_ASSERT(!key.empty());
CAROM_ASSERT(data != 0);
CAROM_ASSERT(nelements > 0);
hsize_t dim[] = { static_cast<hsize_t>(nelements) };
hid_t space = H5Screate_simple(1, dim, 0);
CAROM_ASSERT(space >= 0);
#if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
hid_t dataset = H5Dcreate(d_group_id,
key.c_str(),
H5T_IEEE_F64BE,
space,
H5P_DEFAULT,
H5P_DEFAULT,
H5P_DEFAULT);
#else
hid_t dataset = H5Dcreate(d_group_id,
key.c_str(),
H5T_IEEE_F64BE,
space,
H5P_DEFAULT);
#endif
CAROM_ASSERT(dataset >= 0);
herr_t errf = H5Dwrite(dataset,
H5T_NATIVE_DOUBLE,
H5S_ALL,
H5S_ALL,
H5P_DEFAULT,
data);
CAROM_ASSERT(errf >= 0);
// Write attribute so we know what kind of data this is.
writeAttribute(KEY_DOUBLE_ARRAY, dataset);
errf = H5Sclose(space);
CAROM_ASSERT(errf >= 0);
errf = H5Dclose(dataset);
CAROM_ASSERT(errf >= 0);
#ifndef DEBUG_CHECK_ASSERTIONS
CAROM_NULL_USE(errf);
#endif
}
void
HDFDatabase::getIntegerArray(
const std::string& key,
int* data,
int nelements)
{
CAROM_ASSERT(!key.empty());
#ifndef DEBUG_CHECK_ASSERTIONS
CAROM_NULL_USE(nelements);
#endif
#if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
hid_t dset = H5Dopen(d_group_id, key.c_str(), H5P_DEFAULT);
#else
hid_t dset = H5Dopen(d_group_id, key.c_str());
#endif
CAROM_ASSERT(dset >= 0);
hid_t dspace = H5Dget_space(dset);
CAROM_ASSERT(dspace >= 0);
hsize_t nsel = H5Sget_select_npoints(dspace);
CAROM_ASSERT(static_cast<int>(nsel) == nelements);
herr_t errf;
if (nsel > 0) {
errf = H5Dread(dset, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, data);
CAROM_ASSERT(errf >= 0);
}
errf = H5Sclose(dspace);
CAROM_ASSERT(errf >= 0);
errf = H5Dclose(dset);
CAROM_ASSERT(errf >= 0);
#ifndef DEBUG_CHECK_ASSERTIONS
CAROM_NULL_USE(errf);
#endif
}
void
HDFDatabase::getDoubleArray(
const std::string& key,
double* data,
int nelements)
{
CAROM_ASSERT(!key.empty());
#ifndef DEBUG_CHECK_ASSERTIONS
CAROM_NULL_USE(nelements);
#endif
#if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
hid_t dset = H5Dopen(d_group_id, key.c_str(), H5P_DEFAULT);
#else
hid_t dset = H5Dopen(d_group_id, key.c_str());
#endif
CAROM_ASSERT(dset >= 0);
hid_t dspace = H5Dget_space(dset);
CAROM_ASSERT(dspace >= 0);
hsize_t nsel = H5Sget_select_npoints(dspace);
CAROM_ASSERT(static_cast<int>(nsel) == nelements);
herr_t errf;
if (nsel > 0) {
errf = H5Dread(dset, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, data);
CAROM_ASSERT(errf >= 0);
}
errf = H5Sclose(dspace);
CAROM_ASSERT(errf >= 0);
errf = H5Dclose(dset);
CAROM_ASSERT(errf >= 0);
#ifndef DEBUG_CHECK_ASSERTIONS
CAROM_NULL_USE(errf);
#endif
}
bool
HDFDatabase::isInteger(
const std::string& key)
{
bool is_int = false;
if (!key.empty()) {
#if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
hid_t this_set = H5Dopen(d_group_id, key.c_str(), H5P_DEFAULT);
#else
hid_t this_set = H5Dopen(d_group_id, key.c_str());
#endif
if (this_set > 0) {
int type_key = readAttribute(this_set);
if (type_key == KEY_INT_ARRAY) {
is_int = true;
}
herr_t errf = H5Dclose(this_set);
CAROM_ASSERT(errf >= 0);
#ifndef DEBUG_CHECK_ASSERTIONS
CAROM_NULL_USE(errf);
#endif
}
}
return is_int;
}
bool
HDFDatabase::isDouble(
const std::string& key)
{
bool is_double = false;
if (!key.empty()) {
#if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
hid_t this_set = H5Dopen(d_group_id, key.c_str(), H5P_DEFAULT);
#else
hid_t this_set = H5Dopen(d_group_id, key.c_str());
#endif
if (this_set > 0) {
int type_key = readAttribute(this_set);
if (type_key == KEY_DOUBLE_ARRAY) {
is_double = true;
}
herr_t errf = H5Dclose(this_set);
CAROM_ASSERT(errf >= 0);
#ifndef DEBUG_CHECK_ASSERTIONS
CAROM_NULL_USE(errf);
#endif
}
}
return is_double;
}
void
HDFDatabase::writeAttribute(
int type_key,
hid_t dataset_id)
{
hid_t attr_id = H5Screate(H5S_SCALAR);
CAROM_ASSERT(attr_id >= 0);
#if (H5_VERS_MAJOR > 1) || ((H5_VERS_MAJOR == 1) && (H5_VERS_MINOR > 6))
hid_t attr = H5Acreate(dataset_id,
"Type",
H5T_STD_I8BE,
attr_id,
H5P_DEFAULT,
H5P_DEFAULT);
#else
hid_t attr = H5Acreate(dataset_id,
"Type",
H5T_STD_I8BE,
attr_id,
H5P_DEFAULT);
#endif
CAROM_ASSERT(attr >= 0);
herr_t errf = H5Awrite(attr, H5T_NATIVE_INT, &type_key);
CAROM_ASSERT(errf >= 0);
errf = H5Aclose(attr);
CAROM_ASSERT(errf >= 0);
errf = H5Sclose(attr_id);
CAROM_ASSERT(errf >= 0);
#ifndef DEBUG_CHECK_ASSERTIONS
CAROM_NULL_USE(errf);
#endif
}
int
HDFDatabase::readAttribute(
hid_t dataset_id)
{
hid_t attr = H5Aopen_name(dataset_id, "Type");
CAROM_ASSERT(attr >= 0);
int type_key;
herr_t errf = H5Aread(attr, H5T_NATIVE_INT, &type_key);
CAROM_ASSERT(errf >= 0);
errf = H5Aclose(attr);
CAROM_ASSERT(errf >= 0);
#ifndef DEBUG_CHECK_ASSERTIONS
CAROM_NULL_USE(errf);
#endif
return type_key;
}
}