-
Notifications
You must be signed in to change notification settings - Fork 51
/
Copy pathHDF5IOHandler.cpp
2970 lines (2763 loc) · 99.1 KB
/
HDF5IOHandler.cpp
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
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/* Copyright 2017-2021 Fabian Koller
*
* This file is part of openPMD-api.
*
* openPMD-api is free software: you can redistribute it and/or modify
* it under the terms of of either the GNU General Public License or
* the GNU Lesser General Public License as published by
* the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* openPMD-api is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License and the GNU Lesser General Public License
* for more details.
*
* You should have received a copy of the GNU General Public License
* and the GNU Lesser General Public License along with openPMD-api.
* If not, see <http://www.gnu.org/licenses/>.
*/
#include "openPMD/IO/HDF5/HDF5IOHandler.hpp"
#include "openPMD/IO/HDF5/HDF5IOHandlerImpl.hpp"
#include "openPMD/auxiliary/Environment.hpp"
#include "openPMD/auxiliary/JSON_internal.hpp"
#include "openPMD/auxiliary/Variant.hpp"
#include <optional>
#include <sstream>
#include <stdexcept>
#if openPMD_HAVE_HDF5
#include "openPMD/Datatype.hpp"
#include "openPMD/Error.hpp"
#include "openPMD/IO/HDF5/HDF5Auxiliary.hpp"
#include "openPMD/IO/HDF5/HDF5FilePosition.hpp"
#include "openPMD/IO/IOTask.hpp"
#include "openPMD/auxiliary/Filesystem.hpp"
#include "openPMD/auxiliary/Mpi.hpp"
#include "openPMD/auxiliary/StringManip.hpp"
#include "openPMD/auxiliary/TypeTraits.hpp"
#include "openPMD/backend/Attribute.hpp"
#include <hdf5.h>
#endif
#include <complex>
#include <cstring>
#include <future>
#include <iostream>
#include <stack>
#include <string>
#include <typeinfo>
#include <unordered_map>
#include <utility>
#include <vector>
namespace openPMD
{
#if openPMD_HAVE_HDF5
#if openPMD_USE_VERIFY
#define VERIFY(CONDITION, TEXT) \
{ \
if (!(CONDITION)) \
throw std::runtime_error((TEXT)); \
}
#else
#define VERIFY(CONDITION, TEXT) \
do \
{ \
(void)sizeof(CONDITION); \
} while (0)
#endif
HDF5IOHandlerImpl::HDF5IOHandlerImpl(
AbstractIOHandler *handler,
json::TracingJSON config,
bool do_warn_unused_params)
: AbstractIOHandlerImpl(handler)
, m_datasetTransferProperty{H5P_DEFAULT}
, m_fileAccessProperty{H5P_DEFAULT}
, m_H5T_BOOL_ENUM{H5Tenum_create(H5T_NATIVE_INT8)}
, m_H5T_CFLOAT{H5Tcreate(H5T_COMPOUND, sizeof(float) * 2)}
, m_H5T_CDOUBLE{H5Tcreate(H5T_COMPOUND, sizeof(double) * 2)}
, m_H5T_CLONG_DOUBLE{H5Tcreate(H5T_COMPOUND, sizeof(long double) * 2)}
, m_H5T_LONG_DOUBLE_80_LE{H5Tcopy(H5T_IEEE_F64BE)}
, m_H5T_CLONG_DOUBLE_80_LE{H5Tcreate(H5T_COMPOUND, 16 * 2)}
{
// create a h5py compatible bool type
VERIFY(
m_H5T_BOOL_ENUM >= 0,
"[HDF5] Internal error: Failed to create bool enum");
std::string t{"TRUE"};
std::string f{"FALSE"};
int64_t tVal = 1;
int64_t fVal = 0;
herr_t status;
status = H5Tenum_insert(m_H5T_BOOL_ENUM, t.c_str(), &tVal);
VERIFY(
status == 0, "[HDF5] Internal error: Failed to insert into HDF5 enum");
status = H5Tenum_insert(m_H5T_BOOL_ENUM, f.c_str(), &fVal);
VERIFY(
status == 0, "[HDF5] Internal error: Failed to insert into HDF5 enum");
// create h5py compatible complex types
VERIFY(
m_H5T_CFLOAT >= 0,
"[HDF5] Internal error: Failed to create complex float");
VERIFY(
m_H5T_CDOUBLE >= 0,
"[HDF5] Internal error: Failed to create complex double");
VERIFY(
m_H5T_CLONG_DOUBLE >= 0,
"[HDF5] Internal error: Failed to create complex long double");
H5Tinsert(m_H5T_CFLOAT, "r", 0, H5T_NATIVE_FLOAT);
H5Tinsert(m_H5T_CFLOAT, "i", sizeof(float), H5T_NATIVE_FLOAT);
H5Tinsert(m_H5T_CDOUBLE, "r", 0, H5T_NATIVE_DOUBLE);
H5Tinsert(m_H5T_CDOUBLE, "i", sizeof(double), H5T_NATIVE_DOUBLE);
H5Tinsert(m_H5T_CLONG_DOUBLE, "r", 0, H5T_NATIVE_LDOUBLE);
H5Tinsert(m_H5T_CLONG_DOUBLE, "i", sizeof(long double), H5T_NATIVE_LDOUBLE);
// Create a type that understands 128bit floats with 80 bits of precision
// even on those platforms that do not have it (ARM64, PPC64).
// Otherwise, files created on e.g. AMD64 platforms might not be readable
// on such platforms.
H5Tset_size(m_H5T_LONG_DOUBLE_80_LE, 16);
H5Tset_order(m_H5T_LONG_DOUBLE_80_LE, H5T_ORDER_LE);
H5Tset_precision(m_H5T_LONG_DOUBLE_80_LE, 80);
H5Tset_fields(m_H5T_LONG_DOUBLE_80_LE, 79, 64, 15, 0, 64);
H5Tset_ebias(m_H5T_LONG_DOUBLE_80_LE, 16383);
H5Tset_norm(m_H5T_LONG_DOUBLE_80_LE, H5T_NORM_NONE);
VERIFY(
m_H5T_LONG_DOUBLE_80_LE >= 0,
"[HDF5] Internal error: Failed to create 128-bit long double");
H5Tinsert(m_H5T_CLONG_DOUBLE_80_LE, "r", 0, m_H5T_LONG_DOUBLE_80_LE);
H5Tinsert(m_H5T_CLONG_DOUBLE_80_LE, "i", 16, m_H5T_LONG_DOUBLE_80_LE);
VERIFY(
m_H5T_LONG_DOUBLE_80_LE >= 0,
"[HDF5] Internal error: Failed to create 128-bit complex long double");
// JSON option can overwrite env option:
if (config.json().contains("hdf5"))
{
m_config = config["hdf5"];
{
constexpr char const *const init_json_shadow_str = R"(
{
"dataset": {
"chunks": null
}
})";
auto init_json_shadow = nlohmann::json::parse(init_json_shadow_str);
json::merge(m_config.getShadow(), init_json_shadow);
}
// unused params
if (do_warn_unused_params)
{
auto shadow = m_config.invertShadow();
if (shadow.size() > 0)
{
switch (m_config.originallySpecifiedAs)
{
case json::SupportedLanguages::JSON:
std::cerr
<< "Warning: parts of the backend configuration for "
"HDF5 remain unused:\n"
<< shadow << std::endl;
break;
case json::SupportedLanguages::TOML: {
auto asToml = json::jsonToToml(shadow);
std::cerr
<< "Warning: parts of the backend configuration for "
"HDF5 remain unused:\n"
<< asToml << std::endl;
break;
}
}
}
}
}
#if H5_VERSION_GE(1, 10, 0) && openPMD_HAVE_MPI
auto const hdf5_collective_metadata =
auxiliary::getEnvString("OPENPMD_HDF5_COLLECTIVE_METADATA", "ON");
if (hdf5_collective_metadata == "ON")
m_hdf5_collective_metadata = 1;
else
m_hdf5_collective_metadata = 0;
#endif
}
HDF5IOHandlerImpl::~HDF5IOHandlerImpl()
{
herr_t status;
status = H5Tclose(m_H5T_BOOL_ENUM);
if (status < 0)
std::cerr << "[HDF5] Internal error: Failed to close bool enum\n";
status = H5Tclose(m_H5T_CFLOAT);
if (status < 0)
std::cerr
<< "[HDF5] Internal error: Failed to close complex float type\n";
status = H5Tclose(m_H5T_CDOUBLE);
if (status < 0)
std::cerr
<< "[HDF5] Internal error: Failed to close complex double type\n";
status = H5Tclose(m_H5T_CLONG_DOUBLE);
if (status < 0)
std::cerr << "[HDF5] Internal error: Failed to close complex long "
"double type\n";
status = H5Tclose(m_H5T_LONG_DOUBLE_80_LE);
if (status < 0)
std::cerr
<< "[HDF5] Internal error: Failed to close long double type\n";
status = H5Tclose(m_H5T_CLONG_DOUBLE_80_LE);
if (status < 0)
std::cerr << "[HDF5] Internal error: Failed to close complex long "
"double type\n";
while (!m_openFileIDs.empty())
{
auto file = m_openFileIDs.begin();
status = H5Fclose(*file);
if (status < 0)
std::cerr << "[HDF5] Internal error: Failed to close HDF5 file "
"(serial)\n";
m_openFileIDs.erase(file);
}
if (m_datasetTransferProperty != H5P_DEFAULT)
{
status = H5Pclose(m_datasetTransferProperty);
if (status < 0)
std::cerr << "[HDF5] Internal error: Failed to close HDF5 dataset "
"transfer property\n";
}
if (m_fileAccessProperty != H5P_DEFAULT)
{
status = H5Pclose(m_fileAccessProperty);
if (status < 0)
std::cerr << "[HDF5] Internal error: Failed to close HDF5 file "
"access property\n";
}
}
void HDF5IOHandlerImpl::createFile(
Writable *writable, Parameter<Operation::CREATE_FILE> const ¶meters)
{
if (access::readOnly(m_handler->m_backendAccess))
throw std::runtime_error(
"[HDF5] Creating a file in read-only mode is not possible.");
if (!writable->written)
{
if (!auxiliary::directory_exists(m_handler->directory))
{
bool success = auxiliary::create_directories(m_handler->directory);
VERIFY(
success,
"[HDF5] Internal error: Failed to create directories during "
"HDF5 file creation");
}
std::string name = m_handler->directory + parameters.name;
if (!auxiliary::ends_with(name, ".h5"))
name += ".h5";
unsigned flags{};
switch (m_handler->m_backendAccess)
{
case Access::CREATE:
flags = H5F_ACC_TRUNC;
break;
case Access::APPEND:
if (auxiliary::file_exists(name))
{
flags = H5F_ACC_RDWR;
}
else
{
flags = H5F_ACC_TRUNC;
}
break;
case Access::READ_WRITE:
flags = H5F_ACC_EXCL;
break;
case Access::READ_ONLY:
case Access::READ_LINEAR:
// condition has been checked above
throw std::runtime_error(
"[HDF5] Control flow error in createFile backend access mode.");
}
hid_t id{};
if (flags == H5F_ACC_RDWR)
{
id = H5Fopen(name.c_str(), flags, m_fileAccessProperty);
}
else
{
id = H5Fcreate(
name.c_str(), flags, H5P_DEFAULT, m_fileAccessProperty);
}
VERIFY(id >= 0, "[HDF5] Internal error: Failed to create HDF5 file");
writable->written = true;
writable->abstractFilePosition =
std::make_shared<HDF5FilePosition>("/");
m_fileNames[writable] = name;
m_fileNamesWithID[std::move(name)] = id;
m_openFileIDs.insert(id);
}
}
void HDF5IOHandlerImpl::checkFile(
Writable *, Parameter<Operation::CHECK_FILE> ¶meters)
{
std::string name = m_handler->directory + parameters.name;
if (!auxiliary::ends_with(name, ".h5"))
{
name += ".h5";
}
bool fileExists =
auxiliary::file_exists(name) || auxiliary::directory_exists(name);
#if openPMD_HAVE_MPI
if (m_communicator.has_value())
{
bool fileExistsRes = false;
int status = MPI_Allreduce(
&fileExists,
&fileExistsRes,
1,
MPI_C_BOOL,
MPI_LOR, // logical or
m_communicator.value());
if (status != 0)
{
throw std::runtime_error("MPI Reduction failed!");
}
fileExists = fileExistsRes;
}
#endif
using FileExists = Parameter<Operation::CHECK_FILE>::FileExists;
*parameters.fileExists = fileExists ? FileExists::Yes : FileExists::No;
}
void HDF5IOHandlerImpl::createPath(
Writable *writable, Parameter<Operation::CREATE_PATH> const ¶meters)
{
if (access::readOnly(m_handler->m_backendAccess))
throw std::runtime_error(
"[HDF5] Creating a path in a file opened as read only is not "
"possible.");
hid_t gapl = H5Pcreate(H5P_GROUP_ACCESS);
#if H5_VERSION_GE(1, 10, 0) && openPMD_HAVE_MPI
if (m_hdf5_collective_metadata)
{
H5Pset_all_coll_metadata_ops(gapl, true);
}
#endif
herr_t status;
if (!writable->written)
{
/* Sanitize path */
std::string path = parameters.path;
if (auxiliary::starts_with(path, '/'))
path = auxiliary::replace_first(path, "/", "");
if (!auxiliary::ends_with(path, '/'))
path += '/';
/* Open H5Object to write into */
Writable *position;
if (writable->parent)
position = writable->parent;
else
position = writable; /* root does not have a parent but might still
have to be written */
File file = getFile(position).value();
hid_t node_id =
H5Gopen(file.id, concrete_h5_file_position(position).c_str(), gapl);
VERIFY(
node_id >= 0,
"[HDF5] Internal error: Failed to open HDF5 group during path "
"creation");
/* Create the path in the file */
std::stack<hid_t> groups;
groups.push(node_id);
for (std::string const &folder : auxiliary::split(path, "/", false))
{
// avoid creation of paths that already exist
htri_t const found =
H5Lexists(groups.top(), folder.c_str(), H5P_DEFAULT);
if (found > 0)
continue;
hid_t group_id = H5Gcreate(
groups.top(),
folder.c_str(),
H5P_DEFAULT,
H5P_DEFAULT,
H5P_DEFAULT);
VERIFY(
group_id >= 0,
"[HDF5] Internal error: Failed to create HDF5 group during "
"path creation");
groups.push(group_id);
}
/* Close the groups */
while (!groups.empty())
{
status = H5Gclose(groups.top());
VERIFY(
status == 0,
"[HDF5] Internal error: Failed to close HDF5 group during path "
"creation");
groups.pop();
}
writable->written = true;
writable->abstractFilePosition =
std::make_shared<HDF5FilePosition>(path);
m_fileNames[writable] = file.name;
}
status = H5Pclose(gapl);
VERIFY(
status == 0,
"[HDF5] Internal error: Failed to close HDF5 property during path "
"creation");
}
void HDF5IOHandlerImpl::createDataset(
Writable *writable, Parameter<Operation::CREATE_DATASET> const ¶meters)
{
if (access::readOnly(m_handler->m_backendAccess))
throw std::runtime_error(
"[HDF5] Creating a dataset in a file opened as read only is not "
"possible.");
if (parameters.joinedDimension.has_value())
{
error::throwOperationUnsupportedInBackend(
"ADIOS1", "Joined Arrays currently only supported in ADIOS2");
}
if (!writable->written)
{
/* Sanitize name */
std::string name = parameters.name;
if (auxiliary::starts_with(name, '/'))
name = auxiliary::replace_first(name, "/", "");
if (auxiliary::ends_with(name, '/'))
name = auxiliary::replace_last(name, "/", "");
std::vector<hsize_t> dims;
std::uint64_t num_elements = 1u;
for (auto const &val : parameters.extent)
{
dims.push_back(static_cast<hsize_t>(val));
num_elements *= val;
}
Datatype d = parameters.dtype;
if (d == Datatype::UNDEFINED)
{
// TODO handle unknown dtype
std::cerr << "[HDF5] Datatype::UNDEFINED caught during dataset "
"creation (serial HDF5)"
<< std::endl;
d = Datatype::BOOL;
}
json::TracingJSON config = [&]() {
if (!m_buffered_dataset_config.has_value())
{
// we are only interested in these values from the global config
constexpr char const *const mask_for_global_conf = R"(
{
"dataset": {
"chunks": null
}
})";
m_buffered_dataset_config = m_config.json();
json::filterByTemplate(
*m_buffered_dataset_config,
nlohmann::json::parse(mask_for_global_conf));
}
auto const &buffered_config = *m_buffered_dataset_config;
auto parsed_config = json::parseOptions(
parameters.options, /* considerFiles = */ false);
if (auto hdf5_config_it = parsed_config.config.find("hdf5");
hdf5_config_it != parsed_config.config.end())
{
auto copy = buffered_config;
json::merge(copy, hdf5_config_it.value());
copy = nlohmann::json{{"hdf5", std::move(copy)}};
parsed_config.config = std::move(copy);
}
else
{
parsed_config.config["hdf5"] = buffered_config;
}
return parsed_config;
}();
// general
bool is_resizable_dataset = false;
if (config.json().contains("resizable"))
{
is_resizable_dataset = config["resizable"].json().get<bool>();
}
using chunking_t = std::vector<hsize_t>;
using compute_chunking_t =
std::variant<chunking_t, std::string /* either "none" or "auto"*/>;
bool chunking_config_from_json = false;
auto throw_chunking_error = [&chunking_config_from_json]() {
if (chunking_config_from_json)
{
throw error::BackendConfigSchema(
{"hdf5", "dataset", "chunks"},
R"(Must be "auto", "none", or a an array of integer.)");
}
else
{
throw error::WrongAPIUsage(
"Environment variable OPENPMD_HDF5_CHUNKS accepts values "
"'auto' and 'none'.");
}
};
compute_chunking_t compute_chunking =
auxiliary::getEnvString("OPENPMD_HDF5_CHUNKS", "auto");
// HDF5 specific
if (config.json().contains("hdf5") &&
config["hdf5"].json().contains("dataset"))
{
json::TracingJSON datasetConfig{config["hdf5"]["dataset"]};
if (datasetConfig.json().contains("chunks"))
{
chunking_config_from_json = true;
auto chunks_json = datasetConfig["chunks"];
if (chunks_json.json().is_string())
{
compute_chunking =
json::asLowerCaseStringDynamic(chunks_json.json())
.value();
}
else if (chunks_json.json().is_array())
{
try
{
compute_chunking =
chunks_json.json().get<std::vector<hsize_t>>();
}
catch (nlohmann::json::type_error const &)
{
throw_chunking_error();
}
}
else
{
throw_chunking_error();
}
}
}
std::optional<chunking_t> chunking = std::visit(
auxiliary::overloaded{
[&](chunking_t &&explicitly_specified)
-> std::optional<chunking_t> {
return std::move(explicitly_specified);
},
[&](std::string const &method_name)
-> std::optional<chunking_t> {
if (method_name == "auto")
{
return getOptimalChunkDims(dims, toBytes(d));
}
else if (method_name == "none")
{
return std::nullopt;
}
else
{
throw_chunking_error();
throw std::runtime_error("Unreachable!");
}
}},
std::move(compute_chunking));
parameters.warnUnusedParameters(
config,
"hdf5",
"Warning: parts of the backend configuration for HDF5 dataset '" +
name + "' remain unused:\n");
hid_t gapl = H5Pcreate(H5P_GROUP_ACCESS);
#if H5_VERSION_GE(1, 10, 0) && openPMD_HAVE_MPI
if (m_hdf5_collective_metadata)
{
H5Pset_all_coll_metadata_ops(gapl, true);
}
#endif
/* Open H5Object to write into */
File file{};
if (auto opt = getFile(writable->parent); opt.has_value())
{
file = opt.value();
}
else
{
throw error::Internal(
"[HDF5] CREATE_DATASET task must have a parent with an "
"associated file.");
}
hid_t node_id =
H5Gopen(file.id, concrete_h5_file_position(writable).c_str(), gapl);
VERIFY(
node_id >= 0,
"[HDF5] Internal error: Failed to open HDF5 group during dataset "
"creation");
if (m_handler->m_backendAccess == Access::APPEND)
{
// The dataset might already exist in the file from a previous run
// We delete it, otherwise we could not create it again with
// possibly different parameters.
if (htri_t link_id = H5Lexists(node_id, name.c_str(), H5P_DEFAULT);
link_id > 0)
{
// This only unlinks, but does not delete the dataset
// Deleting the actual dataset physically is now up to HDF5:
// > when removing an object with H5Ldelete, the HDF5 library
// > should be able to detect and recycle the file space when no
// > other reference to the deleted object exists
// https://github.com/openPMD/openPMD-api/pull/1007#discussion_r867223316
herr_t status = H5Ldelete(node_id, name.c_str(), H5P_DEFAULT);
VERIFY(
status == 0,
"[HDF5] Internal error: Failed to delete old dataset '" +
name + "' from group for overwriting.");
}
else if (link_id < 0)
{
throw std::runtime_error(
"[HDF5] Internal error: Failed to check for link existence "
"of '" +
name + "' inside group for overwriting.");
}
// else: link_id == 0: Link does not exist, nothing to do
}
std::vector<hsize_t> max_dims(dims.begin(), dims.end());
if (is_resizable_dataset)
max_dims.assign(dims.size(), H5F_UNLIMITED);
hid_t space = H5Screate_simple(
static_cast<int>(dims.size()), dims.data(), max_dims.data());
VERIFY(
space >= 0,
"[HDF5] Internal error: Failed to create dataspace during dataset "
"creation");
/* enable chunking on the created dataspace */
hid_t datasetCreationProperty = H5Pcreate(H5P_DATASET_CREATE);
H5Pset_fill_time(datasetCreationProperty, H5D_FILL_TIME_NEVER);
if (num_elements != 0u && chunking.has_value())
{
if (chunking->size() != parameters.extent.size())
{
std::string chunking_printed = [&]() {
if (chunking->empty())
{
return std::string("[]");
}
else
{
std::stringstream s;
auto it = chunking->begin();
auto end = chunking->end();
s << '[' << *it++;
for (; it != end; ++it)
{
s << ", " << *it;
}
s << ']';
return s.str();
}
}();
std::cerr << "[HDF5] Chunking for dataset '" << name
<< "' was specified as " << chunking_printed
<< ", but dataset has dimensionality "
<< parameters.extent.size() << ". Will ignore."
<< std::endl;
}
else
{
herr_t status = H5Pset_chunk(
datasetCreationProperty,
chunking->size(),
chunking->data());
VERIFY(
status == 0,
"[HDF5] Internal error: Failed to set chunk size during "
"dataset creation");
}
}
std::string const &compression = ""; // @todo read from JSON
if (!compression.empty())
std::cerr
<< "[HDF5] Compression not yet implemented in HDF5 backend."
<< std::endl;
/*
{
std::vector< std::string > args = auxiliary::split(compression,
":"); std::string const& format = args[0]; if( (format == "zlib" ||
format == "gzip" || format == "deflate")
&& args.size() == 2 )
{
status = H5Pset_deflate(datasetCreationProperty,
std::stoi(args[1])); VERIFY(status == 0, "[HDF5] Internal error: Failed
to set deflate compression during dataset creation"); } else if( format
== "szip" || format == "nbit" || format == "scaleoffset" ) std::cerr <<
"[HDF5] Compression format " << format
<< " not yet implemented. Data will not be
compressed!"
<< std::endl;
else
std::cerr << "[HDF5] Compression format " << format
<< " unknown. Data will not be compressed!"
<< std::endl;
}
*/
GetH5DataType getH5DataType({
{typeid(bool).name(), m_H5T_BOOL_ENUM},
{typeid(std::complex<float>).name(), m_H5T_CFLOAT},
{typeid(std::complex<double>).name(), m_H5T_CDOUBLE},
{typeid(std::complex<long double>).name(), m_H5T_CLONG_DOUBLE},
});
Attribute a(0);
a.dtype = d;
hid_t datatype = getH5DataType(a);
VERIFY(
datatype >= 0,
"[HDF5] Internal error: Failed to get HDF5 datatype during dataset "
"creation");
hid_t group_id = H5Dcreate(
node_id,
name.c_str(),
datatype,
space,
H5P_DEFAULT,
datasetCreationProperty,
H5P_DEFAULT);
VERIFY(
group_id >= 0,
"[HDF5] Internal error: Failed to create HDF5 group during dataset "
"creation");
herr_t status;
status = H5Dclose(group_id);
VERIFY(
status == 0,
"[HDF5] Internal error: Failed to close HDF5 dataset during "
"dataset creation");
status = H5Tclose(datatype);
VERIFY(
status == 0,
"[HDF5] Internal error: Failed to close HDF5 datatype during "
"dataset creation");
status = H5Pclose(datasetCreationProperty);
VERIFY(
status == 0,
"[HDF5] Internal error: Failed to close HDF5 dataset creation "
"property during dataset creation");
status = H5Sclose(space);
VERIFY(
status == 0,
"[HDF5] Internal error: Failed to close HDF5 dataset space during "
"dataset creation");
status = H5Gclose(node_id);
VERIFY(
status == 0,
"[HDF5] Internal error: Failed to close HDF5 group during dataset "
"creation");
status = H5Pclose(gapl);
VERIFY(
status == 0,
"[HDF5] Internal error: Failed to close HDF5 property during "
"dataset creation");
writable->written = true;
writable->abstractFilePosition =
std::make_shared<HDF5FilePosition>(name);
m_fileNames[writable] = file.name;
}
}
void HDF5IOHandlerImpl::extendDataset(
Writable *writable, Parameter<Operation::EXTEND_DATASET> const ¶meters)
{
if (access::readOnly(m_handler->m_backendAccess))
throw std::runtime_error(
"[HDF5] Extending a dataset in a file opened as read only is not "
"possible.");
if (!writable->written)
throw std::runtime_error(
"[HDF5] Extending an unwritten Dataset is not possible.");
auto res = getFile(writable);
if (!res)
res = getFile(writable->parent);
hid_t dataset_id = H5Dopen(
res.value().id,
concrete_h5_file_position(writable).c_str(),
H5P_DEFAULT);
VERIFY(
dataset_id >= 0,
"[HDF5] Internal error: Failed to open HDF5 dataset during dataset "
"extension");
// Datasets may only be extended if they have chunked layout, so let's see
// whether this one does
{
hid_t dataset_space = H5Dget_space(dataset_id);
int ndims = H5Sget_simple_extent_ndims(dataset_space);
VERIFY(
ndims >= 0,
"[HDF5]: Internal error: Failed to retrieve dimensionality of "
"dataset during dataset read.");
hid_t propertyList = H5Dget_create_plist(dataset_id);
std::vector<hsize_t> chunkExtent(ndims, 0);
int chunkDimensionality =
H5Pget_chunk(propertyList, ndims, chunkExtent.data());
if (chunkDimensionality < 0)
{
throw std::runtime_error(
"[HDF5] Cannot extend datasets unless written with chunked "
"layout.");
}
}
std::vector<hsize_t> size;
for (auto const &val : parameters.extent)
size.push_back(static_cast<hsize_t>(val));
herr_t status;
status = H5Dset_extent(dataset_id, size.data());
VERIFY(
status == 0,
"[HDF5] Internal error: Failed to extend HDF5 dataset during dataset "
"extension");
status = H5Dclose(dataset_id);
VERIFY(
status == 0,
"[HDF5] Internal error: Failed to close HDF5 dataset during dataset "
"extension");
}
void HDF5IOHandlerImpl::availableChunks(
Writable *writable, Parameter<Operation::AVAILABLE_CHUNKS> ¶meters)
{
auto fname = m_fileNames.find(writable);
VERIFY(
fname != m_fileNames.end(), "[HDF5] File name not found in writable");
auto fid = m_fileNamesWithID.find(fname->second);
VERIFY(
fid != m_fileNamesWithID.end(),
"[HDF5] File ID not found with file name");
hid_t dataset_id = H5Dopen(
fid->second, concrete_h5_file_position(writable).c_str(), H5P_DEFAULT);
VERIFY(
dataset_id >= 0,
"[HDF5] Internal error: Failed to open HDF5 dataset during dataset "
"read");
hid_t dataset_space = H5Dget_space(dataset_id);
int ndims = H5Sget_simple_extent_ndims(dataset_space);
VERIFY(
ndims >= 0,
"[HDF5]: Internal error: Failed to retrieve dimensionality of "
"dataset "
"during dataset read.");
// // now let's figure out whether this one has chunks
// hid_t propertyList = H5Dget_create_plist( dataset_id );
// std::vector< hsize_t > chunkExtent( ndims, 0 );
// int chunkDimensionality =
// H5Pget_chunk( propertyList, ndims, chunkExtent.data() );
// if( chunkDimensionality >= 0 )
// {
// /*
// * so, the dataset indeed has chunks
// * alas, this backend doesn't write chunks, so for now, reading them
// * is unimplemented
// *
// * https://hdf5.io/develop/group___h5_d.html#gaccff213d3e0765b86f66d08dd9959807
// * May or may not be helpful if implementing this properly one day.
// */
// }
std::vector<hsize_t> dims(ndims, 0);
// return value is equal to ndims
H5Sget_simple_extent_dims(dataset_space, dims.data(), nullptr);
Offset offset(ndims, 0);
Extent extent;
extent.reserve(ndims);
for (auto e : dims)
{
extent.push_back(e);
}
parameters.chunks->push_back(
WrittenChunkInfo(std::move(offset), std::move(extent)));
herr_t status;
status = H5Sclose(dataset_space);
VERIFY(
status == 0,
"[HDF5] Internal error: Failed to close HDF5 dataset space during "
"availableChunks task");
status = H5Dclose(dataset_id);
VERIFY(
status == 0,
"[HDF5] Internal error: Failed to close HDF5 dataset during "
"availableChunks task");
}
void HDF5IOHandlerImpl::openFile(
Writable *writable, Parameter<Operation::OPEN_FILE> ¶meters)
{
if (!auxiliary::directory_exists(m_handler->directory))
throw error::ReadError(
error::AffectedObject::File,
error::Reason::Inaccessible,
"HDF5",
"Supplied directory is not valid: " + m_handler->directory);
std::string name = m_handler->directory + parameters.name;
if (!auxiliary::ends_with(name, ".h5"))
name += ".h5";
// this may (intentionally) overwrite
m_fileNames[writable] = name;
// check if file already open
auto search = m_fileNamesWithID.find(name);
if (search != m_fileNamesWithID.end())
{
return;
}
unsigned flags;
Access at = m_handler->m_backendAccess;
if (access::readOnly(at))
flags = H5F_ACC_RDONLY;
/*
* Within the HDF5 backend, APPEND and READ_WRITE mode are
* equivalent, but the openPMD frontend exposes no reading
* functionality in APPEND mode.
*/
else
flags = H5F_ACC_RDWR;
hid_t file_id;
file_id = H5Fopen(name.c_str(), flags, m_fileAccessProperty);
if (file_id < 0)
throw error::ReadError(
error::AffectedObject::File,
error::Reason::Inaccessible,
"HDF5",
"Failed to open HDF5 file " + name);
writable->written = true;
writable->abstractFilePosition = std::make_shared<HDF5FilePosition>("/");