Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
hdf5.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2019 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_hdf5_h
16#define dealii_hdf5_h
17
18#include <deal.II/base/config.h>
19
20#ifdef DEAL_II_WITH_HDF5
21
23
25
26# include <hdf5.h>
27
28# include <numeric>
29
31
32// It is necessary to turn clang-format off in order to maintain the Doxygen
33// links because they are longer than 80 characters
34// clang-format off
343// clang-format on
344namespace HDF5
345{
350 {
351 protected:
356 HDF5Object(const std::string &name, const bool mpi);
357
358 public:
371 template <typename T>
372 T
373 get_attribute(const std::string &attr_name) const;
374
387 template <typename T>
388 void
389 set_attribute(const std::string &attr_name, const T value);
390
396 std::string
397 get_name() const;
398
399 protected:
405 const std::string name;
406
414 std::shared_ptr<hid_t> hdf5_reference;
415
419 const bool mpi;
420 };
421
425 class DataSet : public HDF5Object
426 {
427 friend class Group;
428
429 protected:
434 DataSet(const std::string &name, const hid_t &parent_group_id, bool mpi);
435
440 DataSet(const std::string &name,
441 const hid_t &parent_group_id,
442 const std::vector<hsize_t> &dimensions,
443 const std::shared_ptr<hid_t> &t_type,
444 const bool mpi);
445
446 public:
464 template <typename Container>
465 Container
466 read();
467
501 template <typename Container>
502 Container
503 read_selection(const std::vector<hsize_t> &coordinates);
504
505 // clang-format off
539 // clang-format on
540 template <typename Container>
541 Container
542 read_hyperslab(const std::vector<hsize_t> &offset,
543 const std::vector<hsize_t> &count);
544
577 template <typename Container>
578 Container
579 read_hyperslab(const std::vector<hsize_t> &data_dimensions,
580 const std::vector<hsize_t> &offset,
581 const std::vector<hsize_t> &stride,
582 const std::vector<hsize_t> &count,
583 const std::vector<hsize_t> &block);
584
596 template <typename number>
597 void
598 read_none();
599
618 template <typename Container>
619 void
620 write(const Container &data);
621
650 template <typename Container>
651 void
652 write_selection(const Container &data,
653 const std::vector<hsize_t> &coordinates);
654
655 // clang-format off
681 // clang-format on
682 template <typename Container>
683 void
684 write_hyperslab(const Container &data,
685 const std::vector<hsize_t> &offset,
686 const std::vector<hsize_t> &count);
687
720 template <typename Container>
721 void
722 write_hyperslab(const Container &data,
723 const std::vector<hsize_t> &data_dimensions,
724 const std::vector<hsize_t> &offset,
725 const std::vector<hsize_t> &stride,
726 const std::vector<hsize_t> &count,
727 const std::vector<hsize_t> &block);
728
746 template <typename number>
747 void
748 write_none();
749
766 bool
767 get_query_io_mode() const;
768
772 void
773 set_query_io_mode(const bool new_query_io_mode);
774
789 std::string
790 get_io_mode();
791
808 H5D_mpio_actual_io_mode_t
810
829 std::string
831
853 std::uint32_t
855
874 std::string
876
897 std::uint32_t
899
905 std::vector<hsize_t>
906 get_dimensions() const;
907
911 unsigned int
912 get_size() const;
913
917 unsigned int
918 get_rank() const;
919
920 private:
924 unsigned int rank;
925
930 std::vector<hsize_t> dimensions;
931
935 std::shared_ptr<hid_t> dataspace;
936
940 unsigned int size;
941
952
956 H5D_mpio_actual_io_mode_t io_mode;
957
964
971 };
972
976 class Group : public HDF5Object
977 {
978 protected:
983 {
987 open,
991 create
992 };
1001 Group(const std::string &name,
1002 const Group &parent_group,
1003 const bool mpi,
1004 const GroupAccessMode mode);
1005
1011 Group(const std::string &name, const bool mpi);
1012
1013 public:
1017 Group
1018 open_group(const std::string &name) const;
1019
1023 Group
1024 create_group(const std::string &name) const;
1025
1029 DataSet
1030 open_dataset(const std::string &name) const;
1031
1042 template <typename number>
1043 DataSet
1044 create_dataset(const std::string &name,
1045 const std::vector<hsize_t> &dimensions) const;
1046
1065 template <typename Container>
1066 void
1067 write_dataset(const std::string &name, const Container &data) const;
1068 };
1069
1073 class File : public Group
1074 {
1075 public:
1080 {
1084 open,
1088 create
1089 };
1090
1096 File(const std::string &name, const FileAccessMode mode);
1097
1105 File(const std::string &name,
1106 const FileAccessMode mode,
1107 const MPI_Comm mpi_communicator);
1108
1109 private:
1117 File(const std::string &name,
1118 const FileAccessMode mode,
1119 const bool mpi,
1120 const MPI_Comm mpi_communicator);
1121 };
1122
1123 namespace internal
1124 {
1136 template <typename number>
1137 std::shared_ptr<hid_t>
1139
1149 template <typename number>
1150 std::vector<hsize_t>
1151 get_container_dimensions(const std::vector<number> &data);
1152
1157 template <typename number>
1158 std::vector<hsize_t>
1160
1165 template <typename number>
1166 std::vector<hsize_t>
1168
1173 template <typename number>
1174 unsigned int
1175 get_container_size(const std::vector<number> &data);
1176
1181 template <typename number>
1182 unsigned int
1184
1189 template <typename number>
1190 unsigned int
1192
1211 template <typename Container>
1212 std::enable_if_t<
1213 std::is_same_v<Container, std::vector<typename Container::value_type>>,
1214 Container>
1215 initialize_container(const std::vector<hsize_t> &dimensions);
1216
1220 template <typename Container>
1221 std::enable_if_t<
1222 std::is_same_v<Container, Vector<typename Container::value_type>>,
1223 Container>
1224 initialize_container(const std::vector<hsize_t> &dimensions);
1225
1229 template <typename Container>
1230 std::enable_if_t<
1231 std::is_same_v<Container, FullMatrix<typename Container::value_type>>,
1232 Container>
1233 initialize_container(const std::vector<hsize_t> &dimensions);
1234
1241 inline void
1242 set_plist(hid_t &plist, const bool mpi);
1243
1252 inline void
1253 release_plist(hid_t &plist,
1254 H5D_mpio_actual_io_mode_t &io_mode,
1255 std::uint32_t &local_no_collective_cause,
1256 std::uint32_t &global_no_collective_cause,
1257 const bool mpi,
1258 const bool query_io_mode);
1259
1263 inline std::string
1264 no_collective_cause_to_string(const std::uint32_t no_collective_cause);
1265 } // namespace internal
1266
1267
1268
1269 // definitions
1270
1271 namespace internal
1272 {
1273 template <typename number>
1274 std::shared_ptr<hid_t>
1276 {
1277 static_assert(std::is_same_v<number, float> ||
1278 std::is_same_v<number, double> ||
1279 std::is_same_v<number, int> ||
1280 std::is_same_v<number, bool> ||
1281 std::is_same_v<number, unsigned int> ||
1282 std::is_same_v<number, std::complex<float>> ||
1283 std::is_same_v<number, std::complex<double>>,
1284 "The data type you are trying to get the HDF5 tag for "
1285 "is not supported by this function.");
1286
1287 if (std::is_same_v<number, float>)
1288 {
1289 return std::make_shared<hid_t>(H5T_NATIVE_FLOAT);
1290 }
1291 else if (std::is_same_v<number, double>)
1292 {
1293 return std::make_shared<hid_t>(H5T_NATIVE_DOUBLE);
1294 }
1295 else if (std::is_same_v<number, int>)
1296 {
1297 return std::make_shared<hid_t>(H5T_NATIVE_INT);
1298 }
1299 else if (std::is_same_v<number, unsigned int>)
1300 {
1301 return std::make_shared<hid_t>(H5T_NATIVE_UINT);
1302 }
1303 else if (std::is_same_v<number, bool>)
1304 {
1305 return std::make_shared<hid_t>(H5T_NATIVE_HBOOL);
1306 }
1307 else if (std::is_same_v<number, std::complex<float>>)
1308 {
1309 std::shared_ptr<hid_t> t_type =
1310 std::shared_ptr<hid_t>(new hid_t, [](hid_t *pointer) {
1311 // Release the HDF5 resource
1312 const herr_t ret = H5Tclose(*pointer);
1313 AssertNothrow(ret >= 0, ExcInternalError());
1314 (void)ret;
1315 delete pointer;
1316 });
1317
1318 *t_type = H5Tcreate(H5T_COMPOUND, sizeof(std::complex<float>));
1319 // The C++ standards committee agreed to mandate that the storage
1320 // format used for the std::complex type be binary-compatible with
1321 // the C99 type, i.e. an array T[2] with consecutive real [0] and
1322 // imaginary [1] parts.
1323 herr_t ret = H5Tinsert(*t_type, "r", 0, H5T_NATIVE_FLOAT);
1324 Assert(ret >= 0, ExcInternalError());
1325 ret = H5Tinsert(*t_type, "i", sizeof(float), H5T_NATIVE_FLOAT);
1326 Assert(ret >= 0, ExcInternalError());
1327 (void)ret;
1328 return t_type;
1329 }
1330 else if (std::is_same_v<number, std::complex<double>>)
1331 {
1332 std::shared_ptr<hid_t> t_type =
1333 std::shared_ptr<hid_t>(new hid_t, [](hid_t *pointer) {
1334 // Release the HDF5 resource
1335 const herr_t ret = H5Tclose(*pointer);
1336 AssertNothrow(ret >= 0, ExcInternalError());
1337 (void)ret;
1338 delete pointer;
1339 });
1340 *t_type = H5Tcreate(H5T_COMPOUND, sizeof(std::complex<double>));
1341 // The C++ standards committee agreed to mandate that the storage
1342 // format used for the std::complex type be binary-compatible with
1343 // the C99 type, i.e. an array T[2] with consecutive real [0] and
1344 // imaginary [1] parts.
1345 herr_t ret = H5Tinsert(*t_type, "r", 0, H5T_NATIVE_DOUBLE);
1346 Assert(ret >= 0, ExcInternalError());
1347 ret = H5Tinsert(*t_type, "i", sizeof(double), H5T_NATIVE_DOUBLE);
1348 Assert(ret >= 0, ExcInternalError());
1349 (void)ret;
1350 return t_type;
1351 }
1352
1353 // The function should not reach this point
1355 return {};
1356 }
1357
1358
1359
1360 template <typename number>
1361 std::vector<hsize_t>
1362 get_container_dimensions(const std::vector<number> &data)
1363 {
1364 std::vector<hsize_t> dimensions = {data.size()};
1365 return dimensions;
1366 }
1367
1368
1369
1370 template <typename number>
1371 std::vector<hsize_t>
1373 {
1374 std::vector<hsize_t> dimensions = {data.size()};
1375 return dimensions;
1376 }
1377
1378
1379
1380 template <typename number>
1381 std::vector<hsize_t>
1383 {
1384 std::vector<hsize_t> dimensions = {data.m(), data.n()};
1385 return dimensions;
1386 }
1387
1388
1389
1390 template <typename number>
1391 unsigned int
1392 get_container_size(const std::vector<number> &data)
1393 {
1394 return static_cast<unsigned int>(data.size());
1395 }
1396
1397
1398
1399 template <typename number>
1400 unsigned int
1402 {
1403 return static_cast<unsigned int>(data.size());
1404 }
1405
1406
1407
1408 template <typename number>
1409 unsigned int
1411 {
1412 return static_cast<unsigned int>(data.m() * data.n());
1413 }
1414
1415
1416
1417 template <typename Container>
1418 std::enable_if_t<
1419 std::is_same_v<Container, std::vector<typename Container::value_type>>,
1420 Container>
1421 initialize_container(const std::vector<hsize_t> &dimensions)
1422 {
1423 return Container(std::accumulate(
1424 dimensions.begin(), dimensions.end(), 1, std::multiplies<int>()));
1425 }
1426
1427
1428
1429 template <typename Container>
1430 std::enable_if_t<
1431 std::is_same_v<Container, Vector<typename Container::value_type>>,
1432 Container>
1433 initialize_container(const std::vector<hsize_t> &dimensions)
1434 {
1435 return Container(std::accumulate(
1436 dimensions.begin(), dimensions.end(), 1, std::multiplies<int>()));
1437 }
1438
1439
1440
1441 template <typename Container>
1442 std::enable_if_t<
1443 std::is_same_v<Container, FullMatrix<typename Container::value_type>>,
1444 Container>
1445 initialize_container(const std::vector<hsize_t> &dimensions)
1446 {
1447 // If the rank is higher than 2, then remove single-dimensional entries
1448 // from the shape defined by dimensions. This is equivalent to the squeeze
1449 // function of python/numpy. For example the following code would convert
1450 // the vector {1,3,1,2} to {3,2}
1451 std::vector<hsize_t> squeezed_dimensions;
1452
1453 if (dimensions.size() > 2)
1454 {
1455 for (const auto &dimension : dimensions)
1456 {
1457 if (dimension > 1)
1458 squeezed_dimensions.push_back(dimension);
1459 }
1460 }
1461 else
1462 {
1463 squeezed_dimensions = dimensions;
1464 }
1465
1466 AssertDimension(squeezed_dimensions.size(), 2);
1467 return Container(squeezed_dimensions[0], squeezed_dimensions[1]);
1468 }
1469
1470
1471 inline void
1472 set_plist(hid_t &plist, const bool mpi)
1473 {
1474 if (mpi)
1475 {
1476# ifdef DEAL_II_WITH_MPI
1477 plist = H5Pcreate(H5P_DATASET_XFER);
1478 Assert(plist >= 0, ExcInternalError());
1479 const herr_t ret = H5Pset_dxpl_mpio(plist, H5FD_MPIO_COLLECTIVE);
1480 (void)ret;
1481 Assert(ret >= 0, ExcInternalError());
1482# else
1484# endif
1485 }
1486 else
1487 {
1488 plist = H5P_DEFAULT;
1489 }
1490
1491 (void)plist;
1492 (void)mpi;
1493 }
1494
1495
1496 inline void
1497 release_plist(hid_t &plist,
1498 H5D_mpio_actual_io_mode_t &io_mode,
1499 std::uint32_t &local_no_collective_cause,
1500 std::uint32_t &global_no_collective_cause,
1501 const bool mpi,
1502 const bool query_io_mode)
1503 {
1504 if (mpi)
1505 {
1506# ifdef DEAL_II_WITH_MPI
1507 herr_t ret;
1508 (void)ret;
1509 if (query_io_mode)
1510 {
1511 ret = H5Pget_mpio_actual_io_mode(plist, &io_mode);
1512 Assert(ret >= 0, ExcInternalError());
1513 ret =
1514 H5Pget_mpio_no_collective_cause(plist,
1515 &local_no_collective_cause,
1516 &global_no_collective_cause);
1517 Assert(ret >= 0, ExcInternalError());
1518 }
1519 ret = H5Pclose(plist);
1520 Assert(ret >= 0, ExcInternalError());
1521# else
1523# endif
1524 }
1525
1526 (void)plist;
1527 (void)io_mode;
1528 (void)local_no_collective_cause;
1529 (void)global_no_collective_cause;
1530 (void)mpi;
1531 (void)query_io_mode;
1532 }
1533
1534
1535 inline std::string
1536 no_collective_cause_to_string(const std::uint32_t no_collective_cause)
1537 {
1538 std::string message;
1539
1540 auto append_to_message = [&message](const char *p) {
1541 if (message.size() > 0)
1542 message += ", ";
1543 message += p;
1544 };
1545
1546 // The first is not a bitmask comparison, the rest are bitmask
1547 // comparisons.
1548 // https://support.hdfgroup.org/HDF5/doc/RM/RM_H5P.html#Property-GetMpioNoCollectiveCause
1549 // See H5Ppublic.h
1550 // Hex codes are used because the HDF5 Group can deprecate some of the
1551 // enum codes. For example the enum code H5D_MPIO_FILTERS is not defined
1552 // in 1.10.2 because it is possible to use compressed datasets with the
1553 // MPI/IO driver.
1554
1555 // H5D_MPIO_COLLECTIVE
1556 if (no_collective_cause == 0x00)
1557 {
1558 append_to_message("H5D_MPIO_COLLECTIVE");
1559 }
1560 // H5D_MPIO_SET_INDEPENDENT
1561 if ((no_collective_cause & 0x01) == 0x01)
1562 {
1563 append_to_message("H5D_MPIO_SET_INDEPENDENT");
1564 }
1565 // H5D_MPIO_DATATYPE_CONVERSION
1566 if ((no_collective_cause & 0x02) == 0x02)
1567 {
1568 append_to_message("H5D_MPIO_DATATYPE_CONVERSION");
1569 }
1570 // H5D_MPIO_DATA_TRANSFORMS
1571 if ((no_collective_cause & 0x04) == 0x04)
1572 {
1573 append_to_message("H5D_MPIO_DATA_TRANSFORMS");
1574 }
1575 // H5D_MPIO_NOT_SIMPLE_OR_SCALAR_DATASPACES
1576 if ((no_collective_cause & 0x10) == 0x10)
1577 {
1578 append_to_message("H5D_MPIO_NOT_SIMPLE_OR_SCALAR_DATASPACES");
1579 }
1580 // H5D_MPIO_NOT_CONTIGUOUS_OR_CHUNKED_DATASET
1581 if ((no_collective_cause & 0x20) == 0x20)
1582 {
1583 append_to_message("H5D_MPIO_NOT_CONTIGUOUS_OR_CHUNKED_DATASET");
1584 }
1585 // H5D_MPIO_FILTERS
1586 if ((no_collective_cause & 0x40) == 0x40)
1587 {
1588 append_to_message("H5D_MPIO_FILTERS");
1589 }
1590 return message;
1591 }
1592 } // namespace internal
1593
1594
1595 template <typename T>
1596 T
1597 HDF5Object::get_attribute(const std::string &attr_name) const
1598 {
1599 const std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<T>();
1600 T value;
1601 hid_t attr;
1602 herr_t ret;
1603
1604 attr = H5Aopen(*hdf5_reference, attr_name.data(), H5P_DEFAULT);
1605 Assert(attr >= 0, ExcMessage("Error at H5Aopen"));
1606 (void)ret;
1607 ret = H5Aread(attr, *t_type, &value);
1608 Assert(ret >= 0, ExcMessage("Error at H5Aread"));
1609 (void)ret;
1610 ret = H5Aclose(attr);
1611 Assert(ret >= 0, ExcMessage("Error at H5Aclose"));
1612
1613 return value;
1614 }
1615
1616
1617
1618 template <>
1619 inline std::string
1620 HDF5Object::get_attribute(const std::string &attr_name) const
1621 {
1622 // Reads a UTF8 variable string
1623 //
1624 // code inspired from
1625 // https://support.hdfgroup.org/ftp/HDF5/examples/misc-examples/vlstratt.c
1626 //
1627 // In the case of a variable length string the user does not have to reserve
1628 // memory for string_out. H5Aread will reserve the memory and the
1629 // user has to free the memory.
1630 //
1631 // Todo:
1632 // - Use H5Dvlen_reclaim instead of free
1633
1634 char *string_out;
1635 hid_t attr;
1636 hid_t type;
1637 herr_t ret;
1638
1639 /* Create a datatype to refer to. */
1640 type = H5Tcopy(H5T_C_S1);
1641 Assert(type >= 0, ExcInternalError());
1642
1643 // Python strings are encoded in UTF8
1644 ret = H5Tset_cset(type, H5T_CSET_UTF8);
1645 Assert(type >= 0, ExcInternalError());
1646
1647 ret = H5Tset_size(type, H5T_VARIABLE);
1648 Assert(ret >= 0, ExcInternalError());
1649
1650 attr = H5Aopen(*hdf5_reference, attr_name.data(), H5P_DEFAULT);
1651 Assert(attr >= 0, ExcInternalError());
1652
1653 ret = H5Aread(attr, type, &string_out);
1654 Assert(ret >= 0, ExcInternalError());
1655
1656 std::string string_value(string_out);
1657 // The memory of the variable length string has to be freed.
1658 // H5Dvlen_reclaim could be also used
1659 free(string_out);
1660 ret = H5Tclose(type);
1661 Assert(ret >= 0, ExcInternalError());
1662
1663 ret = H5Aclose(attr);
1664 Assert(ret >= 0, ExcInternalError());
1665
1666 (void)ret;
1667 return string_value;
1668 }
1669
1670
1671
1672 template <typename T>
1673 void
1674 HDF5Object::set_attribute(const std::string &attr_name, const T value)
1675 {
1676 hid_t attr;
1677 hid_t aid;
1678 herr_t ret;
1679
1680 const std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<T>();
1681
1682
1683 /*
1684 * Create scalar attribute.
1685 */
1686 aid = H5Screate(H5S_SCALAR);
1687 Assert(aid >= 0, ExcMessage("Error at H5Screate"));
1688 attr = H5Acreate2(*hdf5_reference,
1689 attr_name.data(),
1690 *t_type,
1691 aid,
1692 H5P_DEFAULT,
1693 H5P_DEFAULT);
1694 Assert(attr >= 0, ExcMessage("Error at H5Acreate2"));
1695
1696 /*
1697 * Write scalar attribute.
1698 */
1699 ret = H5Awrite(attr, *t_type, &value);
1700 Assert(ret >= 0, ExcMessage("Error at H5Awrite"));
1701
1702 ret = H5Sclose(aid);
1703 Assert(ret >= 0, ExcMessage("Error at H5Sclose"));
1704 ret = H5Aclose(attr);
1705 Assert(ret >= 0, ExcMessage("Error at H5Aclose"));
1706
1707 (void)ret;
1708 }
1709
1710
1711
1712 template <>
1713 inline void
1714 HDF5Object::set_attribute(const std::string &attr_name,
1715 const std::string value) // NOLINT
1716 {
1717 // Writes a UTF8 variable string
1718 //
1719 // code inspired from
1720 // https://support.hdfgroup.org/ftp/HDF5/examples/misc-examples/vlstratt.c
1721
1722 hid_t attr;
1723 hid_t aid;
1724 hid_t t_type;
1725 herr_t ret;
1726
1727 /* Create a datatype to refer to. */
1728 t_type = H5Tcopy(H5T_C_S1);
1729 Assert(t_type >= 0, ExcInternalError());
1730
1731 // Python strings are encoded in UTF8
1732 ret = H5Tset_cset(t_type, H5T_CSET_UTF8);
1733 Assert(t_type >= 0, ExcInternalError());
1734
1735 ret = H5Tset_size(t_type, H5T_VARIABLE);
1736 Assert(ret >= 0, ExcInternalError());
1737
1738 /*
1739 * Create scalar attribute.
1740 */
1741 aid = H5Screate(H5S_SCALAR);
1742 Assert(aid >= 0, ExcMessage("Error at H5Screate"));
1743 attr = H5Acreate2(
1744 *hdf5_reference, attr_name.data(), t_type, aid, H5P_DEFAULT, H5P_DEFAULT);
1745 Assert(attr >= 0, ExcMessage("Error at H5Acreate2"));
1746
1747 /*
1748 * Write scalar attribute.
1749 * In most of the cases H5Awrite and H5Dwrite take a pointer to the data.
1750 * But in the particular case of a variable length string, H5Awrite takes
1751 * the address of the pointer of the string.
1752 */
1753 const char *c_string_value = value.c_str();
1754 ret = H5Awrite(attr, t_type, &c_string_value);
1755 Assert(ret >= 0, ExcInternalError());
1756
1757 ret = H5Sclose(aid);
1758 Assert(ret >= 0, ExcMessage("Error at H5Sclose"));
1759 ret = H5Aclose(attr);
1760 Assert(ret >= 0, ExcMessage("Error at H5Aclose"));
1761
1762 (void)ret;
1763 }
1764
1765
1766
1767 template <typename Container>
1768 Container
1770 {
1771 const std::shared_ptr<hid_t> t_type =
1772 internal::get_hdf5_datatype<typename Container::value_type>();
1773 hid_t plist;
1774 herr_t ret;
1775
1776 Container data = internal::initialize_container<Container>(dimensions);
1777
1778 internal::set_plist(plist, mpi);
1779
1780 ret = H5Dread(*hdf5_reference,
1781 *t_type,
1782 H5S_ALL,
1783 H5S_ALL,
1784 plist,
1785 make_array_view(data).data());
1786 Assert(ret >= 0, ExcInternalError());
1787
1788
1790 io_mode,
1793 mpi,
1795
1796 (void)ret;
1797 return data;
1798 }
1799
1800
1801
1802 template <typename Container>
1803 Container
1804 DataSet::read_selection(const std::vector<hsize_t> &coordinates)
1805 {
1806 Assert(coordinates.size() % rank == 0,
1807 ExcMessage(
1808 "The dimension of coordinates has to be divisible by the rank"));
1809 const std::shared_ptr<hid_t> t_type =
1810 internal::get_hdf5_datatype<typename Container::value_type>();
1811 hid_t plist;
1812 hid_t memory_dataspace;
1813 herr_t ret;
1814
1815 std::vector<hsize_t> data_dimensions{
1816 static_cast<hsize_t>(coordinates.size() / rank)};
1817
1818 Container data = internal::initialize_container<Container>(data_dimensions);
1819
1820 memory_dataspace = H5Screate_simple(1, data_dimensions.data(), nullptr);
1821 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
1822 ret = H5Sselect_elements(*dataspace,
1823 H5S_SELECT_SET,
1824 data.size(),
1825 coordinates.data());
1826 Assert(ret >= 0, ExcMessage("Error at H5Sselect_elements"));
1827
1828 internal::set_plist(plist, mpi);
1829
1830 ret = H5Dread(*hdf5_reference,
1831 *t_type,
1832 memory_dataspace,
1833 *dataspace,
1834 plist,
1835 make_array_view(data).data());
1836 Assert(ret >= 0, ExcMessage("Error at H5Dread"));
1837
1839 io_mode,
1842 mpi,
1844
1845 ret = H5Sclose(memory_dataspace);
1846 Assert(ret >= 0, ExcMessage("Error at H5SClose"));
1847
1848 (void)ret;
1849 return data;
1850 }
1851
1852
1853
1854 template <typename Container>
1855 Container
1856 DataSet::read_hyperslab(const std::vector<hsize_t> &offset,
1857 const std::vector<hsize_t> &count)
1858 {
1859 const std::shared_ptr<hid_t> t_type =
1860 internal::get_hdf5_datatype<typename Container::value_type>();
1861 hid_t plist;
1862 hid_t memory_dataspace;
1863 herr_t ret;
1864
1865 // In this particular overload of read_hyperslab the data_dimensions are
1866 // the same as count
1867 std::vector<hsize_t> data_dimensions = count;
1868
1869 Container data = internal::initialize_container<Container>(data_dimensions);
1870
1871 memory_dataspace =
1872 H5Screate_simple(data_dimensions.size(), data_dimensions.data(), nullptr);
1873 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
1874 ret = H5Sselect_hyperslab(*dataspace,
1875 H5S_SELECT_SET,
1876 offset.data(),
1877 nullptr,
1878 count.data(),
1879 nullptr);
1880 Assert(ret >= 0, ExcMessage("Error at H5Sselect_hyperslab"));
1881
1882 internal::set_plist(plist, mpi);
1883
1884 ret = H5Dread(*hdf5_reference,
1885 *t_type,
1886 memory_dataspace,
1887 *dataspace,
1888 plist,
1889 make_array_view(data).data());
1890 Assert(ret >= 0, ExcMessage("Error at H5Dread"));
1891
1893 io_mode,
1896 mpi,
1898
1899 ret = H5Sclose(memory_dataspace);
1900 Assert(ret >= 0, ExcMessage("Error at H5SClose"));
1901
1902 (void)ret;
1903 return data;
1904 }
1905
1906
1907
1908 template <typename Container>
1909 Container
1910 DataSet::read_hyperslab(const std::vector<hsize_t> &data_dimensions,
1911 const std::vector<hsize_t> &offset,
1912 const std::vector<hsize_t> &stride,
1913 const std::vector<hsize_t> &count,
1914 const std::vector<hsize_t> &block)
1915 {
1916 const std::shared_ptr<hid_t> t_type =
1917 internal::get_hdf5_datatype<typename Container::value_type>();
1918 hid_t plist;
1919 hid_t memory_dataspace;
1920 herr_t ret;
1921
1922 Container data = internal::initialize_container<Container>(data_dimensions);
1923
1924 memory_dataspace =
1925 H5Screate_simple(data_dimensions.size(), data_dimensions.data(), nullptr);
1926 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
1927 ret = H5Sselect_hyperslab(*dataspace,
1928 H5S_SELECT_SET,
1929 offset.data(),
1930 stride.data(),
1931 count.data(),
1932 block.data());
1933 Assert(ret >= 0, ExcMessage("Error at H5Sselect_hyperslab"));
1934
1935 internal::set_plist(plist, mpi);
1936
1937 ret = H5Dread(*hdf5_reference,
1938 *t_type,
1939 memory_dataspace,
1940 *dataspace,
1941 plist,
1942 make_array_view(data).data());
1943 Assert(ret >= 0, ExcMessage("Error at H5Dread"));
1944
1946 io_mode,
1949 mpi,
1951
1952 ret = H5Sclose(memory_dataspace);
1953 Assert(ret >= 0, ExcMessage("Error at H5SClose"));
1954
1955 (void)ret;
1956 return data;
1957 }
1958
1959
1960
1961 template <typename number>
1962 void
1964 {
1965 const std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<number>();
1966 const std::vector<hsize_t> data_dimensions = {0};
1967
1968 hid_t memory_dataspace;
1969 hid_t plist;
1970 herr_t ret;
1971
1972 memory_dataspace = H5Screate_simple(1, data_dimensions.data(), nullptr);
1973 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
1974 ret = H5Sselect_none(*dataspace);
1975 Assert(ret >= 0, ExcMessage("H5Sselect_none"));
1976
1977 internal::set_plist(plist, mpi);
1978
1979 // The pointer of data can safely be nullptr, see the discussion at the HDF5
1980 // forum:
1981 // https://forum.hdfgroup.org/t/parallel-i-o-does-not-support-filters-yet/884/17
1982 ret = H5Dread(
1983 *hdf5_reference, *t_type, memory_dataspace, *dataspace, plist, nullptr);
1984 Assert(ret >= 0, ExcMessage("Error at H5Dread"));
1985
1987 io_mode,
1990 mpi,
1992
1993 ret = H5Sclose(memory_dataspace);
1994 Assert(ret >= 0, ExcMessage("Error at H5SClose"));
1995
1996 (void)ret;
1997 }
1998
1999
2000
2001 template <typename Container>
2002 void
2003 DataSet::write(const Container &data)
2004 {
2006 const std::shared_ptr<hid_t> t_type =
2007 internal::get_hdf5_datatype<typename Container::value_type>();
2008 hid_t plist;
2009 herr_t ret;
2010
2011 internal::set_plist(plist, mpi);
2012
2013 ret = H5Dwrite(*hdf5_reference,
2014 *t_type,
2015 H5S_ALL,
2016 H5S_ALL,
2017 plist,
2018 make_array_view(data).data());
2019 Assert(ret >= 0, ExcMessage("Error at H5Dwrite"));
2020
2022 io_mode,
2025 mpi,
2027
2028 (void)ret;
2029 }
2030
2031
2032
2033 template <typename Container>
2034 void
2035 DataSet::write_selection(const Container &data,
2036 const std::vector<hsize_t> &coordinates)
2037 {
2038 AssertDimension(coordinates.size(), data.size() * rank);
2039 const std::shared_ptr<hid_t> t_type =
2040 internal::get_hdf5_datatype<typename Container::value_type>();
2041 const std::vector<hsize_t> data_dimensions =
2043
2044 hid_t memory_dataspace;
2045 hid_t plist;
2046 herr_t ret;
2047
2048
2049 memory_dataspace = H5Screate_simple(1, data_dimensions.data(), nullptr);
2050 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
2051 ret = H5Sselect_elements(*dataspace,
2052 H5S_SELECT_SET,
2053 data.size(),
2054 coordinates.data());
2055 Assert(ret >= 0, ExcMessage("Error at H5Sselect_elements"));
2056
2057 internal::set_plist(plist, mpi);
2058
2059 ret = H5Dwrite(*hdf5_reference,
2060 *t_type,
2061 memory_dataspace,
2062 *dataspace,
2063 plist,
2064 make_array_view(data).data());
2065 Assert(ret >= 0, ExcMessage("Error at H5Dwrite"));
2066
2068 io_mode,
2071 mpi,
2073
2074 ret = H5Sclose(memory_dataspace);
2075 Assert(ret >= 0, ExcMessage("Error at H5SClose"));
2076
2077 (void)ret;
2078 }
2079
2080
2081
2082 template <typename Container>
2083 void
2084 DataSet::write_hyperslab(const Container &data,
2085 const std::vector<hsize_t> &offset,
2086 const std::vector<hsize_t> &count)
2087 {
2088 AssertDimension(std::accumulate(count.begin(),
2089 count.end(),
2090 1,
2091 std::multiplies<unsigned int>()),
2093 const std::shared_ptr<hid_t> t_type =
2094 internal::get_hdf5_datatype<typename Container::value_type>();
2095 // In this particular overload of write_hyperslab the data_dimensions are
2096 // the same as count
2097 const std::vector<hsize_t> &data_dimensions = count;
2098
2099 hid_t memory_dataspace;
2100 hid_t plist;
2101 herr_t ret;
2102
2103 memory_dataspace =
2104 H5Screate_simple(data_dimensions.size(), data_dimensions.data(), nullptr);
2105 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
2106 ret = H5Sselect_hyperslab(*dataspace,
2107 H5S_SELECT_SET,
2108 offset.data(),
2109 nullptr,
2110 count.data(),
2111 nullptr);
2112 Assert(ret >= 0, ExcMessage("Error at H5Sselect_hyperslab"));
2113
2114 internal::set_plist(plist, mpi);
2115
2116 ret = H5Dwrite(*hdf5_reference,
2117 *t_type,
2118 memory_dataspace,
2119 *dataspace,
2120 plist,
2121 make_array_view(data).data());
2122 Assert(ret >= 0, ExcMessage("Error at H5Dwrite"));
2123
2125 io_mode,
2128 mpi,
2130
2131 ret = H5Sclose(memory_dataspace);
2132 Assert(ret >= 0, ExcMessage("Error at H5Sclose"));
2133
2134 (void)ret;
2135 }
2136
2137
2138
2139 template <typename Container>
2140 void
2141 DataSet::write_hyperslab(const Container &data,
2142 const std::vector<hsize_t> &data_dimensions,
2143 const std::vector<hsize_t> &offset,
2144 const std::vector<hsize_t> &stride,
2145 const std::vector<hsize_t> &count,
2146 const std::vector<hsize_t> &block)
2147 {
2148 const std::shared_ptr<hid_t> t_type =
2149 internal::get_hdf5_datatype<typename Container::value_type>();
2150
2151 hid_t memory_dataspace;
2152 hid_t plist;
2153 herr_t ret;
2154
2155 memory_dataspace =
2156 H5Screate_simple(data_dimensions.size(), data_dimensions.data(), nullptr);
2157 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
2158 ret = H5Sselect_hyperslab(*dataspace,
2159 H5S_SELECT_SET,
2160 offset.data(),
2161 stride.data(),
2162 count.data(),
2163 block.data());
2164 Assert(ret >= 0, ExcMessage("Error at H5Sselect_hyperslab"));
2165
2166 internal::set_plist(plist, mpi);
2167
2168 ret = H5Dwrite(*hdf5_reference,
2169 *t_type,
2170 memory_dataspace,
2171 *dataspace,
2172 plist,
2173 make_array_view(data).data());
2174 Assert(ret >= 0, ExcMessage("Error at H5Dwrite"));
2175
2177 io_mode,
2180 mpi,
2182
2183 ret = H5Sclose(memory_dataspace);
2184 Assert(ret >= 0, ExcMessage("Error at H5Sclose"));
2185
2186 (void)ret;
2187 }
2188
2189
2190
2191 template <typename number>
2192 void
2194 {
2195 std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<number>();
2196 std::vector<hsize_t> data_dimensions = {0};
2197
2198 hid_t memory_dataspace;
2199 hid_t plist;
2200 herr_t ret;
2201
2202 memory_dataspace = H5Screate_simple(1, data_dimensions.data(), nullptr);
2203 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
2204 ret = H5Sselect_none(*dataspace);
2205 Assert(ret >= 0, ExcMessage("Error at H5PSselect_none"));
2206
2207 internal::set_plist(plist, mpi);
2208
2209 // The pointer of data can safely be nullptr, see the discussion at the HDF5
2210 // forum:
2211 // https://forum.hdfgroup.org/t/parallel-i-o-does-not-support-filters-yet/884/17
2212 ret = H5Dwrite(
2213 *hdf5_reference, *t_type, memory_dataspace, *dataspace, plist, nullptr);
2214 Assert(ret >= 0, ExcMessage("Error at H5Dwrite"));
2215
2217 io_mode,
2220 mpi,
2222
2223 ret = H5Sclose(memory_dataspace);
2224 Assert(ret >= 0, ExcMessage("Error at H5Sclose"));
2225
2226 (void)ret;
2227 }
2228
2229
2230
2231 template <typename number>
2232 DataSet
2233 Group::create_dataset(const std::string &name,
2234 const std::vector<hsize_t> &dimensions) const
2235 {
2236 std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<number>();
2237 return {name, *hdf5_reference, dimensions, t_type, mpi};
2238 }
2239
2240
2241
2242 template <typename Container>
2243 void
2244 Group::write_dataset(const std::string &name, const Container &data) const
2245 {
2246 std::vector<hsize_t> dimensions = internal::get_container_dimensions(data);
2247 auto dataset =
2248 create_dataset<typename Container::value_type>(name, dimensions);
2249 dataset.write(data);
2250 }
2251} // namespace HDF5
2252
2254
2255
2256#endif // DEAL_II_WITH_HDF5
2257
2258#endif // dealii_hdf5_h
ArrayView< typename std::remove_reference< typename std::iterator_traits< Iterator >::reference >::type, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
Definition array_view.h:949
size_type n() const
size_type m() const
unsigned int rank
Definition hdf5.h:924
unsigned int get_rank() const
Definition hdf5.cc:321
void write(const Container &data)
Definition hdf5.h:2003
bool query_io_mode
Definition hdf5.h:951
std::vector< hsize_t > get_dimensions() const
Definition hdf5.cc:205
std::vector< hsize_t > dimensions
Definition hdf5.h:930
std::uint32_t get_local_no_collective_cause_as_hdf5_type()
Definition hdf5.cc:270
void write_hyperslab(const Container &data, const std::vector< hsize_t > &offset, const std::vector< hsize_t > &count)
Definition hdf5.h:2084
unsigned int get_size() const
Definition hdf5.cc:313
void write_selection(const Container &data, const std::vector< hsize_t > &coordinates)
Definition hdf5.h:2035
void read_none()
Definition hdf5.h:1963
std::shared_ptr< hid_t > dataspace
Definition hdf5.h:935
std::string get_local_no_collective_cause()
Definition hdf5.cc:258
std::string get_io_mode()
Definition hdf5.cc:213
std::uint32_t local_no_collective_cause
Definition hdf5.h:963
std::string get_global_no_collective_cause()
Definition hdf5.cc:282
void write_none()
Definition hdf5.h:2193
Container read()
Definition hdf5.h:1769
bool get_query_io_mode() const
Definition hdf5.cc:305
std::uint32_t global_no_collective_cause
Definition hdf5.h:970
std::uint32_t get_global_no_collective_cause_as_hdf5_type()
Definition hdf5.cc:294
Container read_selection(const std::vector< hsize_t > &coordinates)
Definition hdf5.h:1804
H5D_mpio_actual_io_mode_t get_io_mode_as_hdf5_type()
Definition hdf5.cc:247
void set_query_io_mode(const bool new_query_io_mode)
Definition hdf5.cc:197
Container read_hyperslab(const std::vector< hsize_t > &offset, const std::vector< hsize_t > &count)
Definition hdf5.h:1856
unsigned int size
Definition hdf5.h:940
H5D_mpio_actual_io_mode_t io_mode
Definition hdf5.h:956
FileAccessMode
Definition hdf5.h:1080
DataSet open_dataset(const std::string &name) const
Definition hdf5.cc:388
Group open_group(const std::string &name) const
Definition hdf5.cc:372
DataSet create_dataset(const std::string &name, const std::vector< hsize_t > &dimensions) const
Definition hdf5.h:2233
Group create_group(const std::string &name) const
Definition hdf5.cc:380
GroupAccessMode
Definition hdf5.h:983
void write_dataset(const std::string &name, const Container &data) const
Definition hdf5.h:2244
std::shared_ptr< hid_t > hdf5_reference
Definition hdf5.h:414
std::string get_name() const
Definition hdf5.cc:81
T get_attribute(const std::string &attr_name) const
Definition hdf5.h:1597
const std::string name
Definition hdf5.h:405
void set_attribute(const std::string &attr_name, const T value)
Definition hdf5.h:1674
const bool mpi
Definition hdf5.h:419
virtual size_type size() const override
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertNothrow(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
#define DEAL_II_ASSERT_UNREACHABLE()
void release_plist(hid_t &plist, H5D_mpio_actual_io_mode_t &io_mode, std::uint32_t &local_no_collective_cause, std::uint32_t &global_no_collective_cause, const bool mpi, const bool query_io_mode)
Definition hdf5.h:1497
std::string no_collective_cause_to_string(const std::uint32_t no_collective_cause)
Definition hdf5.h:1536
std::enable_if_t< std::is_same_v< Container, std::vector< typename Container::value_type > >, Container > initialize_container(const std::vector< hsize_t > &dimensions)
Definition hdf5.h:1421
unsigned int get_container_size(const std::vector< number > &data)
Definition hdf5.h:1392
void set_plist(hid_t &plist, const bool mpi)
Definition hdf5.h:1472
std::shared_ptr< hid_t > get_hdf5_datatype()
Definition hdf5.h:1275
std::vector< hsize_t > get_container_dimensions(const std::vector< number > &data)
Definition hdf5.h:1362
Definition hdf5.h:345