deal.II version 9.7.0
\(\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 - 2025 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
342// clang-format on
343namespace HDF5
344{
349 {
350 protected:
355 HDF5Object(const std::string &name, const bool mpi);
356
357 public:
370 template <typename T>
371 T
372 get_attribute(const std::string &attr_name) const;
373
386 template <typename T>
387 void
388 set_attribute(const std::string &attr_name, const T value);
389
395 std::string
396 get_name() const;
397
398 protected:
404 const std::string name;
405
413 std::shared_ptr<hid_t> hdf5_reference;
414
418 const bool mpi;
419 };
420
424 class DataSet : public HDF5Object
425 {
426 friend class Group;
427
428 protected:
433 DataSet(const std::string &name, const hid_t &parent_group_id, bool mpi);
434
439 DataSet(const std::string &name,
440 const hid_t &parent_group_id,
441 const std::vector<hsize_t> &dimensions,
442 const std::shared_ptr<hid_t> &t_type,
443 const bool mpi);
444
445 public:
463 template <typename Container>
464 Container
465 read();
466
500 template <typename Container>
501 Container
502 read_selection(const std::vector<hsize_t> &coordinates);
503
504 // clang-format off
538 // clang-format on
539 template <typename Container>
540 Container
541 read_hyperslab(const std::vector<hsize_t> &offset,
542 const std::vector<hsize_t> &count);
543
576 template <typename Container>
577 Container
578 read_hyperslab(const std::vector<hsize_t> &data_dimensions,
579 const std::vector<hsize_t> &offset,
580 const std::vector<hsize_t> &stride,
581 const std::vector<hsize_t> &count,
582 const std::vector<hsize_t> &block);
583
595 template <typename number>
596 void
597 read_none();
598
617 template <typename Container>
618 void
619 write(const Container &data);
620
649 template <typename Container>
650 void
651 write_selection(const Container &data,
652 const std::vector<hsize_t> &coordinates);
653
654 // clang-format off
680 // clang-format on
681 template <typename Container>
682 void
683 write_hyperslab(const Container &data,
684 const std::vector<hsize_t> &offset,
685 const std::vector<hsize_t> &count);
686
719 template <typename Container>
720 void
721 write_hyperslab(const Container &data,
722 const std::vector<hsize_t> &data_dimensions,
723 const std::vector<hsize_t> &offset,
724 const std::vector<hsize_t> &stride,
725 const std::vector<hsize_t> &count,
726 const std::vector<hsize_t> &block);
727
745 template <typename number>
746 void
747 write_none();
748
765 bool
766 get_query_io_mode() const;
767
771 void
772 set_query_io_mode(const bool new_query_io_mode);
773
788 std::string
789 get_io_mode();
790
807 H5D_mpio_actual_io_mode_t
809
828 std::string
830
852 std::uint32_t
854
873 std::string
875
896 std::uint32_t
898
904 std::vector<hsize_t>
905 get_dimensions() const;
906
910 unsigned int
911 get_size() const;
912
916 unsigned int
917 get_rank() const;
918
919 private:
923 unsigned int rank;
924
929 std::vector<hsize_t> dimensions;
930
934 std::shared_ptr<hid_t> dataspace;
935
939 unsigned int size;
940
951
955 H5D_mpio_actual_io_mode_t io_mode;
956
963
970 };
971
975 class Group : public HDF5Object
976 {
977 protected:
982 {
991 };
992
1000 Group(const std::string &name,
1001 const Group &parent_group,
1002 const bool mpi,
1003 const GroupAccessMode mode);
1004
1010 Group(const std::string &name, const bool mpi);
1011
1012 public:
1016 Group
1017 open_group(const std::string &name) const;
1018
1022 Group
1023 create_group(const std::string &name) const;
1024
1028 DataSet
1029 open_dataset(const std::string &name) const;
1030
1041 template <typename number>
1042 DataSet
1043 create_dataset(const std::string &name,
1044 const std::vector<hsize_t> &dimensions) const;
1045
1064 template <typename Container>
1065 void
1066 write_dataset(const std::string &name, const Container &data) const;
1067 };
1068
1072 class File : public Group
1073 {
1074 public:
1079 {
1088 };
1089
1095 File(const std::string &name, const FileAccessMode mode);
1096
1104 File(const std::string &name,
1105 const FileAccessMode mode,
1106 const MPI_Comm mpi_communicator);
1107
1108 private:
1116 File(const std::string &name,
1117 const FileAccessMode mode,
1118 const bool mpi,
1119 const MPI_Comm mpi_communicator);
1120 };
1121
1122 namespace internal
1123 {
1135 template <typename number>
1136 std::shared_ptr<hid_t>
1138
1148 template <typename number>
1149 std::vector<hsize_t>
1150 get_container_dimensions(const std::vector<number> &data);
1151
1156 template <typename number>
1157 std::vector<hsize_t>
1159
1164 template <typename number>
1165 std::vector<hsize_t>
1167
1172 template <typename number>
1173 unsigned int
1174 get_container_size(const std::vector<number> &data);
1175
1180 template <typename number>
1181 unsigned int
1183
1188 template <typename number>
1189 unsigned int
1191
1210 template <typename Container>
1211 std::enable_if_t<
1212 std::is_same_v<Container, std::vector<typename Container::value_type>>,
1213 Container>
1214 initialize_container(const std::vector<hsize_t> &dimensions);
1215
1219 template <typename Container>
1220 std::enable_if_t<
1221 std::is_same_v<Container, Vector<typename Container::value_type>>,
1222 Container>
1223 initialize_container(const std::vector<hsize_t> &dimensions);
1224
1228 template <typename Container>
1229 std::enable_if_t<
1230 std::is_same_v<Container, FullMatrix<typename Container::value_type>>,
1231 Container>
1232 initialize_container(const std::vector<hsize_t> &dimensions);
1233
1240 inline void
1241 set_plist(hid_t &plist, const bool mpi);
1242
1251 inline void
1252 release_plist(hid_t &plist,
1253 H5D_mpio_actual_io_mode_t &io_mode,
1254 std::uint32_t &local_no_collective_cause,
1255 std::uint32_t &global_no_collective_cause,
1256 const bool mpi,
1257 const bool query_io_mode);
1258
1262 inline std::string
1263 no_collective_cause_to_string(const std::uint32_t no_collective_cause);
1264 } // namespace internal
1265
1266
1267
1268 // definitions
1269
1270 namespace internal
1271 {
1272 template <typename number>
1273 std::shared_ptr<hid_t>
1275 {
1276 static_assert(std::is_same_v<number, float> ||
1277 std::is_same_v<number, double> ||
1278 std::is_same_v<number, int> ||
1279 std::is_same_v<number, bool> ||
1280 std::is_same_v<number, unsigned int> ||
1281 std::is_same_v<number, std::complex<float>> ||
1282 std::is_same_v<number, std::complex<double>>,
1283 "The data type you are trying to get the HDF5 tag for "
1284 "is not supported by this function.");
1285
1286 if (std::is_same_v<number, float>)
1287 {
1288 return std::make_shared<hid_t>(H5T_NATIVE_FLOAT);
1289 }
1290 else if (std::is_same_v<number, double>)
1291 {
1292 return std::make_shared<hid_t>(H5T_NATIVE_DOUBLE);
1293 }
1294 else if (std::is_same_v<number, int>)
1295 {
1296 return std::make_shared<hid_t>(H5T_NATIVE_INT);
1297 }
1298 else if (std::is_same_v<number, unsigned int>)
1299 {
1300 return std::make_shared<hid_t>(H5T_NATIVE_UINT);
1301 }
1302 else if (std::is_same_v<number, bool>)
1303 {
1304 return std::make_shared<hid_t>(H5T_NATIVE_HBOOL);
1305 }
1306 else if (std::is_same_v<number, std::complex<float>>)
1307 {
1308 std::shared_ptr<hid_t> t_type =
1309 std::shared_ptr<hid_t>(new hid_t, [](hid_t *pointer) {
1310 // Release the HDF5 resource
1311 const herr_t ret = H5Tclose(*pointer);
1312 AssertNothrow(ret >= 0, ExcInternalError());
1313 delete pointer;
1314 });
1315
1316 *t_type = H5Tcreate(H5T_COMPOUND, sizeof(std::complex<float>));
1317 // The C++ standards committee agreed to mandate that the storage
1318 // format used for the std::complex type be binary-compatible with
1319 // the C99 type, i.e. an array T[2] with consecutive real [0] and
1320 // imaginary [1] parts.
1321 herr_t ret = H5Tinsert(*t_type, "r", 0, H5T_NATIVE_FLOAT);
1322 Assert(ret >= 0, ExcInternalError());
1323 ret = H5Tinsert(*t_type, "i", sizeof(float), H5T_NATIVE_FLOAT);
1324 Assert(ret >= 0, ExcInternalError());
1325 (void)ret;
1326 return t_type;
1327 }
1328 else if (std::is_same_v<number, std::complex<double>>)
1329 {
1330 std::shared_ptr<hid_t> t_type =
1331 std::shared_ptr<hid_t>(new hid_t, [](hid_t *pointer) {
1332 // Release the HDF5 resource
1333 const herr_t ret = H5Tclose(*pointer);
1334 AssertNothrow(ret >= 0, ExcInternalError());
1335 delete pointer;
1336 });
1337 *t_type = H5Tcreate(H5T_COMPOUND, sizeof(std::complex<double>));
1338 // The C++ standards committee agreed to mandate that the storage
1339 // format used for the std::complex type be binary-compatible with
1340 // the C99 type, i.e. an array T[2] with consecutive real [0] and
1341 // imaginary [1] parts.
1342 herr_t ret = H5Tinsert(*t_type, "r", 0, H5T_NATIVE_DOUBLE);
1343 Assert(ret >= 0, ExcInternalError());
1344 ret = H5Tinsert(*t_type, "i", sizeof(double), H5T_NATIVE_DOUBLE);
1345 Assert(ret >= 0, ExcInternalError());
1346 (void)ret;
1347 return t_type;
1348 }
1349
1350 // The function should not reach this point
1352 return {};
1353 }
1354
1355
1356
1357 template <typename number>
1358 std::vector<hsize_t>
1359 get_container_dimensions(const std::vector<number> &data)
1360 {
1361 std::vector<hsize_t> dimensions = {data.size()};
1362 return dimensions;
1363 }
1364
1365
1366
1367 template <typename number>
1368 std::vector<hsize_t>
1370 {
1371 std::vector<hsize_t> dimensions = {data.size()};
1372 return dimensions;
1373 }
1374
1375
1376
1377 template <typename number>
1378 std::vector<hsize_t>
1380 {
1381 std::vector<hsize_t> dimensions = {data.m(), data.n()};
1382 return dimensions;
1383 }
1384
1385
1386
1387 template <typename number>
1388 unsigned int
1389 get_container_size(const std::vector<number> &data)
1390 {
1391 return static_cast<unsigned int>(data.size());
1392 }
1393
1394
1395
1396 template <typename number>
1397 unsigned int
1399 {
1400 return static_cast<unsigned int>(data.size());
1401 }
1402
1403
1404
1405 template <typename number>
1406 unsigned int
1408 {
1409 return static_cast<unsigned int>(data.m() * data.n());
1410 }
1411
1412
1413
1414 template <typename Container>
1415 std::enable_if_t<
1416 std::is_same_v<Container, std::vector<typename Container::value_type>>,
1417 Container>
1418 initialize_container(const std::vector<hsize_t> &dimensions)
1419 {
1420 return Container(std::accumulate(
1421 dimensions.begin(), dimensions.end(), 1, std::multiplies<int>()));
1422 }
1423
1424
1425
1426 template <typename Container>
1427 std::enable_if_t<
1428 std::is_same_v<Container, Vector<typename Container::value_type>>,
1429 Container>
1430 initialize_container(const std::vector<hsize_t> &dimensions)
1431 {
1432 return Container(std::accumulate(
1433 dimensions.begin(), dimensions.end(), 1, std::multiplies<int>()));
1434 }
1435
1436
1437
1438 template <typename Container>
1439 std::enable_if_t<
1440 std::is_same_v<Container, FullMatrix<typename Container::value_type>>,
1441 Container>
1442 initialize_container(const std::vector<hsize_t> &dimensions)
1443 {
1444 // If the rank is higher than 2, then remove single-dimensional entries
1445 // from the shape defined by dimensions. This is equivalent to the squeeze
1446 // function of python/numpy. For example the following code would convert
1447 // the vector {1,3,1,2} to {3,2}
1448 std::vector<hsize_t> squeezed_dimensions;
1449
1450 if (dimensions.size() > 2)
1451 {
1452 for (const auto &dimension : dimensions)
1453 {
1454 if (dimension > 1)
1455 squeezed_dimensions.push_back(dimension);
1456 }
1457 }
1458 else
1459 {
1460 squeezed_dimensions = dimensions;
1461 }
1462
1463 AssertDimension(squeezed_dimensions.size(), 2);
1464 return Container(squeezed_dimensions[0], squeezed_dimensions[1]);
1465 }
1466
1467
1468 inline void
1469 set_plist(hid_t &plist, const bool mpi)
1470 {
1471 if (mpi)
1472 {
1473# ifdef DEAL_II_WITH_MPI
1474 plist = H5Pcreate(H5P_DATASET_XFER);
1475 Assert(plist >= 0, ExcInternalError());
1476 const herr_t ret = H5Pset_dxpl_mpio(plist, H5FD_MPIO_COLLECTIVE);
1477 (void)ret;
1478 Assert(ret >= 0, ExcInternalError());
1479# else
1481# endif
1482 }
1483 else
1484 {
1485 plist = H5P_DEFAULT;
1486 }
1487
1488 (void)plist;
1489 (void)mpi;
1490 }
1491
1492
1493 inline void
1494 release_plist(hid_t &plist,
1495 H5D_mpio_actual_io_mode_t &io_mode,
1496 std::uint32_t &local_no_collective_cause,
1497 std::uint32_t &global_no_collective_cause,
1498 const bool mpi,
1499 const bool query_io_mode)
1500 {
1501 if (mpi)
1502 {
1503# ifdef DEAL_II_WITH_MPI
1504 herr_t ret;
1505 (void)ret;
1506 if (query_io_mode)
1507 {
1508 ret = H5Pget_mpio_actual_io_mode(plist, &io_mode);
1509 Assert(ret >= 0, ExcInternalError());
1510 ret =
1511 H5Pget_mpio_no_collective_cause(plist,
1512 &local_no_collective_cause,
1513 &global_no_collective_cause);
1514 Assert(ret >= 0, ExcInternalError());
1515 }
1516 ret = H5Pclose(plist);
1517 Assert(ret >= 0, ExcInternalError());
1518# else
1520# endif
1521 }
1522
1523 (void)plist;
1524 (void)io_mode;
1525 (void)local_no_collective_cause;
1526 (void)global_no_collective_cause;
1527 (void)mpi;
1528 (void)query_io_mode;
1529 }
1530
1531
1532 inline std::string
1533 no_collective_cause_to_string(const std::uint32_t no_collective_cause)
1534 {
1535 std::string message;
1536
1537 auto append_to_message = [&message](const char *p) {
1538 if (message.size() > 0)
1539 message += ", ";
1540 message += p;
1541 };
1542
1543 // The first is not a bitmask comparison, the rest are bitmask
1544 // comparisons.
1545 // https://support.hdfgroup.org/HDF5/doc/RM/RM_H5P.html#Property-GetMpioNoCollectiveCause
1546 // See H5Ppublic.h
1547 // Hex codes are used because the HDF5 Group can deprecate some of the
1548 // enum codes. For example the enum code H5D_MPIO_FILTERS is not defined
1549 // in 1.10.2 because it is possible to use compressed datasets with the
1550 // MPI/IO driver.
1551
1552 // H5D_MPIO_COLLECTIVE
1553 if (no_collective_cause == 0x00)
1554 {
1555 append_to_message("H5D_MPIO_COLLECTIVE");
1556 }
1557 // H5D_MPIO_SET_INDEPENDENT
1558 if ((no_collective_cause & 0x01) == 0x01)
1559 {
1560 append_to_message("H5D_MPIO_SET_INDEPENDENT");
1561 }
1562 // H5D_MPIO_DATATYPE_CONVERSION
1563 if ((no_collective_cause & 0x02) == 0x02)
1564 {
1565 append_to_message("H5D_MPIO_DATATYPE_CONVERSION");
1566 }
1567 // H5D_MPIO_DATA_TRANSFORMS
1568 if ((no_collective_cause & 0x04) == 0x04)
1569 {
1570 append_to_message("H5D_MPIO_DATA_TRANSFORMS");
1571 }
1572 // H5D_MPIO_NOT_SIMPLE_OR_SCALAR_DATASPACES
1573 if ((no_collective_cause & 0x10) == 0x10)
1574 {
1575 append_to_message("H5D_MPIO_NOT_SIMPLE_OR_SCALAR_DATASPACES");
1576 }
1577 // H5D_MPIO_NOT_CONTIGUOUS_OR_CHUNKED_DATASET
1578 if ((no_collective_cause & 0x20) == 0x20)
1579 {
1580 append_to_message("H5D_MPIO_NOT_CONTIGUOUS_OR_CHUNKED_DATASET");
1581 }
1582 // H5D_MPIO_FILTERS
1583 if ((no_collective_cause & 0x40) == 0x40)
1584 {
1585 append_to_message("H5D_MPIO_FILTERS");
1586 }
1587 return message;
1588 }
1589 } // namespace internal
1590
1591
1592 template <typename T>
1593 T
1594 HDF5Object::get_attribute(const std::string &attr_name) const
1595 {
1596 const std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<T>();
1597 T value;
1598 hid_t attr;
1599 herr_t ret;
1600
1601 attr = H5Aopen(*hdf5_reference, attr_name.data(), H5P_DEFAULT);
1602 Assert(attr >= 0, ExcMessage("Error at H5Aopen"));
1603 (void)ret;
1604 ret = H5Aread(attr, *t_type, &value);
1605 Assert(ret >= 0, ExcMessage("Error at H5Aread"));
1606 (void)ret;
1607 ret = H5Aclose(attr);
1608 Assert(ret >= 0, ExcMessage("Error at H5Aclose"));
1609
1610 return value;
1611 }
1612
1613
1614
1615 template <>
1616 inline std::string
1617 HDF5Object::get_attribute(const std::string &attr_name) const
1618 {
1619 // Reads a UTF8 variable string
1620 //
1621 // code inspired from
1622 // https://support.hdfgroup.org/ftp/HDF5/examples/misc-examples/vlstratt.c
1623 //
1624 // In the case of a variable length string the user does not have to reserve
1625 // memory for string_out. H5Aread will reserve the memory and the
1626 // user has to free the memory.
1627 //
1628 // Todo:
1629 // - Use H5Dvlen_reclaim instead of free
1630
1631 char *string_out;
1632 hid_t attr;
1633 hid_t type;
1634 herr_t ret;
1635
1636 /* Create a datatype to refer to. */
1637 type = H5Tcopy(H5T_C_S1);
1638 Assert(type >= 0, ExcInternalError());
1639
1640 // Python strings are encoded in UTF8
1641 ret = H5Tset_cset(type, H5T_CSET_UTF8);
1642 Assert(type >= 0, ExcInternalError());
1643
1644 ret = H5Tset_size(type, H5T_VARIABLE);
1645 Assert(ret >= 0, ExcInternalError());
1646
1647 attr = H5Aopen(*hdf5_reference, attr_name.data(), H5P_DEFAULT);
1648 Assert(attr >= 0, ExcInternalError());
1649
1650 ret = H5Aread(attr, type, &string_out);
1651 Assert(ret >= 0, ExcInternalError());
1652
1653 std::string string_value(string_out);
1654 // The memory of the variable length string has to be freed.
1655 // H5Dvlen_reclaim could be also used
1656 std::free(string_out);
1657 ret = H5Tclose(type);
1658 Assert(ret >= 0, ExcInternalError());
1659
1660 ret = H5Aclose(attr);
1661 Assert(ret >= 0, ExcInternalError());
1662
1663 (void)ret;
1664 return string_value;
1665 }
1666
1667
1668
1669 template <typename T>
1670 void
1671 HDF5Object::set_attribute(const std::string &attr_name, const T value)
1672 {
1673 hid_t attr;
1674 hid_t aid;
1675 herr_t ret;
1676
1677 const std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<T>();
1678
1679
1680 /*
1681 * Create scalar attribute.
1682 */
1683 aid = H5Screate(H5S_SCALAR);
1684 Assert(aid >= 0, ExcMessage("Error at H5Screate"));
1685 attr = H5Acreate2(*hdf5_reference,
1686 attr_name.data(),
1687 *t_type,
1688 aid,
1689 H5P_DEFAULT,
1690 H5P_DEFAULT);
1691 Assert(attr >= 0, ExcMessage("Error at H5Acreate2"));
1692
1693 /*
1694 * Write scalar attribute.
1695 */
1696 ret = H5Awrite(attr, *t_type, &value);
1697 Assert(ret >= 0, ExcMessage("Error at H5Awrite"));
1698
1699 ret = H5Sclose(aid);
1700 Assert(ret >= 0, ExcMessage("Error at H5Sclose"));
1701 ret = H5Aclose(attr);
1702 Assert(ret >= 0, ExcMessage("Error at H5Aclose"));
1703
1704 (void)ret;
1705 }
1706
1707
1708
1709 template <>
1710 inline void
1711 HDF5Object::set_attribute(const std::string &attr_name,
1712 const std::string value) // NOLINT
1713 {
1714 // Writes a UTF8 variable string
1715 //
1716 // code inspired from
1717 // https://support.hdfgroup.org/ftp/HDF5/examples/misc-examples/vlstratt.c
1718
1719 hid_t attr;
1720 hid_t aid;
1721 hid_t t_type;
1722 herr_t ret;
1723
1724 /* Create a datatype to refer to. */
1725 t_type = H5Tcopy(H5T_C_S1);
1726 Assert(t_type >= 0, ExcInternalError());
1727
1728 // Python strings are encoded in UTF8
1729 ret = H5Tset_cset(t_type, H5T_CSET_UTF8);
1730 Assert(t_type >= 0, ExcInternalError());
1731
1732 ret = H5Tset_size(t_type, H5T_VARIABLE);
1733 Assert(ret >= 0, ExcInternalError());
1734
1735 /*
1736 * Create scalar attribute.
1737 */
1738 aid = H5Screate(H5S_SCALAR);
1739 Assert(aid >= 0, ExcMessage("Error at H5Screate"));
1740 attr = H5Acreate2(
1741 *hdf5_reference, attr_name.data(), t_type, aid, H5P_DEFAULT, H5P_DEFAULT);
1742 Assert(attr >= 0, ExcMessage("Error at H5Acreate2"));
1743
1744 /*
1745 * Write scalar attribute.
1746 * In most of the cases H5Awrite and H5Dwrite take a pointer to the data.
1747 * But in the particular case of a variable length string, H5Awrite takes
1748 * the address of the pointer of the string.
1749 */
1750 const char *c_string_value = value.c_str();
1751 ret = H5Awrite(attr, t_type, &c_string_value);
1752 Assert(ret >= 0, ExcInternalError());
1753
1754 ret = H5Sclose(aid);
1755 Assert(ret >= 0, ExcMessage("Error at H5Sclose"));
1756 ret = H5Aclose(attr);
1757 Assert(ret >= 0, ExcMessage("Error at H5Aclose"));
1758
1759 (void)ret;
1760 }
1761
1762
1763
1764 template <typename Container>
1765 Container
1767 {
1768 const std::shared_ptr<hid_t> t_type =
1770 hid_t plist;
1771 herr_t ret;
1772
1774
1775 internal::set_plist(plist, mpi);
1776
1777 ret = H5Dread(*hdf5_reference,
1778 *t_type,
1779 H5S_ALL,
1780 H5S_ALL,
1781 plist,
1782 make_array_view(data).data());
1783 Assert(ret >= 0, ExcInternalError());
1784
1785
1787 io_mode,
1790 mpi,
1792
1793 (void)ret;
1794 return data;
1795 }
1796
1797
1798
1799 template <typename Container>
1800 Container
1801 DataSet::read_selection(const std::vector<hsize_t> &coordinates)
1802 {
1803 Assert(coordinates.size() % rank == 0,
1804 ExcMessage(
1805 "The dimension of coordinates has to be divisible by the rank"));
1806 const std::shared_ptr<hid_t> t_type =
1808 hid_t plist;
1809 hid_t memory_dataspace;
1810 herr_t ret;
1811
1812 std::vector<hsize_t> data_dimensions{
1813 static_cast<hsize_t>(coordinates.size() / rank)};
1814
1815 Container data = internal::initialize_container<Container>(data_dimensions);
1816
1817 memory_dataspace = H5Screate_simple(1, data_dimensions.data(), nullptr);
1818 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
1819 ret = H5Sselect_elements(*dataspace,
1820 H5S_SELECT_SET,
1821 data.size(),
1822 coordinates.data());
1823 Assert(ret >= 0, ExcMessage("Error at H5Sselect_elements"));
1824
1825 internal::set_plist(plist, mpi);
1826
1827 ret = H5Dread(*hdf5_reference,
1828 *t_type,
1829 memory_dataspace,
1830 *dataspace,
1831 plist,
1832 make_array_view(data).data());
1833 Assert(ret >= 0, ExcMessage("Error at H5Dread"));
1834
1836 io_mode,
1839 mpi,
1841
1842 ret = H5Sclose(memory_dataspace);
1843 Assert(ret >= 0, ExcMessage("Error at H5SClose"));
1844
1845 (void)ret;
1846 return data;
1847 }
1848
1849
1850
1851 template <typename Container>
1852 Container
1853 DataSet::read_hyperslab(const std::vector<hsize_t> &offset,
1854 const std::vector<hsize_t> &count)
1855 {
1856 const std::shared_ptr<hid_t> t_type =
1858 hid_t plist;
1859 hid_t memory_dataspace;
1860 herr_t ret;
1861
1862 // In this particular overload of read_hyperslab the data_dimensions are
1863 // the same as count
1864 std::vector<hsize_t> data_dimensions = count;
1865
1866 Container data = internal::initialize_container<Container>(data_dimensions);
1867
1868 memory_dataspace =
1869 H5Screate_simple(data_dimensions.size(), data_dimensions.data(), nullptr);
1870 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
1871 ret = H5Sselect_hyperslab(*dataspace,
1872 H5S_SELECT_SET,
1873 offset.data(),
1874 nullptr,
1875 count.data(),
1876 nullptr);
1877 Assert(ret >= 0, ExcMessage("Error at H5Sselect_hyperslab"));
1878
1879 internal::set_plist(plist, mpi);
1880
1881 ret = H5Dread(*hdf5_reference,
1882 *t_type,
1883 memory_dataspace,
1884 *dataspace,
1885 plist,
1886 make_array_view(data).data());
1887 Assert(ret >= 0, ExcMessage("Error at H5Dread"));
1888
1890 io_mode,
1893 mpi,
1895
1896 ret = H5Sclose(memory_dataspace);
1897 Assert(ret >= 0, ExcMessage("Error at H5SClose"));
1898
1899 (void)ret;
1900 return data;
1901 }
1902
1903
1904
1905 template <typename Container>
1906 Container
1907 DataSet::read_hyperslab(const std::vector<hsize_t> &data_dimensions,
1908 const std::vector<hsize_t> &offset,
1909 const std::vector<hsize_t> &stride,
1910 const std::vector<hsize_t> &count,
1911 const std::vector<hsize_t> &block)
1912 {
1913 const std::shared_ptr<hid_t> t_type =
1915 hid_t plist;
1916 hid_t memory_dataspace;
1917 herr_t ret;
1918
1919 Container data = internal::initialize_container<Container>(data_dimensions);
1920
1921 memory_dataspace =
1922 H5Screate_simple(data_dimensions.size(), data_dimensions.data(), nullptr);
1923 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
1924 ret = H5Sselect_hyperslab(*dataspace,
1925 H5S_SELECT_SET,
1926 offset.data(),
1927 stride.data(),
1928 count.data(),
1929 block.data());
1930 Assert(ret >= 0, ExcMessage("Error at H5Sselect_hyperslab"));
1931
1932 internal::set_plist(plist, mpi);
1933
1934 ret = H5Dread(*hdf5_reference,
1935 *t_type,
1936 memory_dataspace,
1937 *dataspace,
1938 plist,
1939 make_array_view(data).data());
1940 Assert(ret >= 0, ExcMessage("Error at H5Dread"));
1941
1943 io_mode,
1946 mpi,
1948
1949 ret = H5Sclose(memory_dataspace);
1950 Assert(ret >= 0, ExcMessage("Error at H5SClose"));
1951
1952 (void)ret;
1953 return data;
1954 }
1955
1956
1957
1958 template <typename number>
1959 void
1961 {
1962 const std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<number>();
1963 const std::vector<hsize_t> data_dimensions = {0};
1964
1965 hid_t memory_dataspace;
1966 hid_t plist;
1967 herr_t ret;
1968
1969 memory_dataspace = H5Screate_simple(1, data_dimensions.data(), nullptr);
1970 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
1971 ret = H5Sselect_none(*dataspace);
1972 Assert(ret >= 0, ExcMessage("H5Sselect_none"));
1973
1974 internal::set_plist(plist, mpi);
1975
1976 // The pointer of data can safely be nullptr, see the discussion at the HDF5
1977 // forum:
1978 // https://forum.hdfgroup.org/t/parallel-i-o-does-not-support-filters-yet/884/17
1979 ret = H5Dread(
1980 *hdf5_reference, *t_type, memory_dataspace, *dataspace, plist, nullptr);
1981 Assert(ret >= 0, ExcMessage("Error at H5Dread"));
1982
1984 io_mode,
1987 mpi,
1989
1990 ret = H5Sclose(memory_dataspace);
1991 Assert(ret >= 0, ExcMessage("Error at H5SClose"));
1992
1993 (void)ret;
1994 }
1995
1996
1997
1998 template <typename Container>
1999 void
2000 DataSet::write(const Container &data)
2001 {
2003 const std::shared_ptr<hid_t> t_type =
2005 hid_t plist;
2006 herr_t ret;
2007
2008 internal::set_plist(plist, mpi);
2009
2010 ret = H5Dwrite(*hdf5_reference,
2011 *t_type,
2012 H5S_ALL,
2013 H5S_ALL,
2014 plist,
2015 make_array_view(data).data());
2016 Assert(ret >= 0, ExcMessage("Error at H5Dwrite"));
2017
2019 io_mode,
2022 mpi,
2024
2025 (void)ret;
2026 }
2027
2028
2029
2030 template <typename Container>
2031 void
2032 DataSet::write_selection(const Container &data,
2033 const std::vector<hsize_t> &coordinates)
2034 {
2035 AssertDimension(coordinates.size(), data.size() * rank);
2036 const std::shared_ptr<hid_t> t_type =
2038 const std::vector<hsize_t> data_dimensions =
2040
2041 hid_t memory_dataspace;
2042 hid_t plist;
2043 herr_t ret;
2044
2045
2046 memory_dataspace = H5Screate_simple(1, data_dimensions.data(), nullptr);
2047 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
2048 ret = H5Sselect_elements(*dataspace,
2049 H5S_SELECT_SET,
2050 data.size(),
2051 coordinates.data());
2052 Assert(ret >= 0, ExcMessage("Error at H5Sselect_elements"));
2053
2054 internal::set_plist(plist, mpi);
2055
2056 ret = H5Dwrite(*hdf5_reference,
2057 *t_type,
2058 memory_dataspace,
2059 *dataspace,
2060 plist,
2061 make_array_view(data).data());
2062 Assert(ret >= 0, ExcMessage("Error at H5Dwrite"));
2063
2065 io_mode,
2068 mpi,
2070
2071 ret = H5Sclose(memory_dataspace);
2072 Assert(ret >= 0, ExcMessage("Error at H5SClose"));
2073
2074 (void)ret;
2075 }
2076
2077
2078
2079 template <typename Container>
2080 void
2081 DataSet::write_hyperslab(const Container &data,
2082 const std::vector<hsize_t> &offset,
2083 const std::vector<hsize_t> &count)
2084 {
2085 AssertDimension(std::accumulate(count.begin(),
2086 count.end(),
2087 1,
2088 std::multiplies<unsigned int>()),
2090 const std::shared_ptr<hid_t> t_type =
2092 // In this particular overload of write_hyperslab the data_dimensions are
2093 // the same as count
2094 const std::vector<hsize_t> &data_dimensions = count;
2095
2096 hid_t memory_dataspace;
2097 hid_t plist;
2098 herr_t ret;
2099
2100 memory_dataspace =
2101 H5Screate_simple(data_dimensions.size(), data_dimensions.data(), nullptr);
2102 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
2103 ret = H5Sselect_hyperslab(*dataspace,
2104 H5S_SELECT_SET,
2105 offset.data(),
2106 nullptr,
2107 count.data(),
2108 nullptr);
2109 Assert(ret >= 0, ExcMessage("Error at H5Sselect_hyperslab"));
2110
2111 internal::set_plist(plist, mpi);
2112
2113 ret = H5Dwrite(*hdf5_reference,
2114 *t_type,
2115 memory_dataspace,
2116 *dataspace,
2117 plist,
2118 make_array_view(data).data());
2119 Assert(ret >= 0, ExcMessage("Error at H5Dwrite"));
2120
2122 io_mode,
2125 mpi,
2127
2128 ret = H5Sclose(memory_dataspace);
2129 Assert(ret >= 0, ExcMessage("Error at H5Sclose"));
2130
2131 (void)ret;
2132 }
2133
2134
2135
2136 template <typename Container>
2137 void
2138 DataSet::write_hyperslab(const Container &data,
2139 const std::vector<hsize_t> &data_dimensions,
2140 const std::vector<hsize_t> &offset,
2141 const std::vector<hsize_t> &stride,
2142 const std::vector<hsize_t> &count,
2143 const std::vector<hsize_t> &block)
2144 {
2145 const std::shared_ptr<hid_t> t_type =
2147
2148 hid_t memory_dataspace;
2149 hid_t plist;
2150 herr_t ret;
2151
2152 memory_dataspace =
2153 H5Screate_simple(data_dimensions.size(), data_dimensions.data(), nullptr);
2154 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
2155 ret = H5Sselect_hyperslab(*dataspace,
2156 H5S_SELECT_SET,
2157 offset.data(),
2158 stride.data(),
2159 count.data(),
2160 block.data());
2161 Assert(ret >= 0, ExcMessage("Error at H5Sselect_hyperslab"));
2162
2163 internal::set_plist(plist, mpi);
2164
2165 ret = H5Dwrite(*hdf5_reference,
2166 *t_type,
2167 memory_dataspace,
2168 *dataspace,
2169 plist,
2170 make_array_view(data).data());
2171 Assert(ret >= 0, ExcMessage("Error at H5Dwrite"));
2172
2174 io_mode,
2177 mpi,
2179
2180 ret = H5Sclose(memory_dataspace);
2181 Assert(ret >= 0, ExcMessage("Error at H5Sclose"));
2182
2183 (void)ret;
2184 }
2185
2186
2187
2188 template <typename number>
2189 void
2191 {
2192 std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<number>();
2193 std::vector<hsize_t> data_dimensions = {0};
2194
2195 hid_t memory_dataspace;
2196 hid_t plist;
2197 herr_t ret;
2198
2199 memory_dataspace = H5Screate_simple(1, data_dimensions.data(), nullptr);
2200 Assert(memory_dataspace >= 0, ExcMessage("Error at H5Screate_simple"));
2201 ret = H5Sselect_none(*dataspace);
2202 Assert(ret >= 0, ExcMessage("Error at H5PSselect_none"));
2203
2204 internal::set_plist(plist, mpi);
2205
2206 // The pointer of data can safely be nullptr, see the discussion at the HDF5
2207 // forum:
2208 // https://forum.hdfgroup.org/t/parallel-i-o-does-not-support-filters-yet/884/17
2209 ret = H5Dwrite(
2210 *hdf5_reference, *t_type, memory_dataspace, *dataspace, plist, nullptr);
2211 Assert(ret >= 0, ExcMessage("Error at H5Dwrite"));
2212
2214 io_mode,
2217 mpi,
2219
2220 ret = H5Sclose(memory_dataspace);
2221 Assert(ret >= 0, ExcMessage("Error at H5Sclose"));
2222
2223 (void)ret;
2224 }
2225
2226
2227
2228 template <typename number>
2229 DataSet
2230 Group::create_dataset(const std::string &name,
2231 const std::vector<hsize_t> &dimensions) const
2232 {
2233 std::shared_ptr<hid_t> t_type = internal::get_hdf5_datatype<number>();
2234 return {name, *hdf5_reference, dimensions, t_type, mpi};
2235 }
2236
2237
2238
2239 template <typename Container>
2240 void
2241 Group::write_dataset(const std::string &name, const Container &data) const
2242 {
2243 std::vector<hsize_t> dimensions = internal::get_container_dimensions(data);
2244 auto dataset =
2246 dataset.write(data);
2247 }
2248} // namespace HDF5
2249
2251
2252
2253#else
2254
2255// Make sure the scripts that create the C++20 module input files have
2256// something to latch on if the preprocessor #ifdef above would
2257// otherwise lead to an empty content of the file.
2260
2261#endif // DEAL_II_WITH_HDF5
2262
2263#endif // dealii_hdf5_h
ArrayView< std::remove_reference_t< typename std::iterator_traits< Iterator >::reference >, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
Definition array_view.h:954
size_type n() const
size_type m() const
unsigned int rank
Definition hdf5.h:923
unsigned int get_rank() const
Definition hdf5.cc:312
void write(const Container &data)
Definition hdf5.h:2000
bool query_io_mode
Definition hdf5.h:950
friend class Group
Definition hdf5.h:426
std::vector< hsize_t > get_dimensions() const
Definition hdf5.cc:196
std::vector< hsize_t > dimensions
Definition hdf5.h:929
std::uint32_t get_local_no_collective_cause_as_hdf5_type()
Definition hdf5.cc:261
void write_hyperslab(const Container &data, const std::vector< hsize_t > &offset, const std::vector< hsize_t > &count)
Definition hdf5.h:2081
unsigned int get_size() const
Definition hdf5.cc:304
void write_selection(const Container &data, const std::vector< hsize_t > &coordinates)
Definition hdf5.h:2032
void read_none()
Definition hdf5.h:1960
std::shared_ptr< hid_t > dataspace
Definition hdf5.h:934
DataSet(const std::string &name, const hid_t &parent_group_id, bool mpi)
Definition hdf5.cc:83
std::string get_local_no_collective_cause()
Definition hdf5.cc:249
std::string get_io_mode()
Definition hdf5.cc:204
std::uint32_t local_no_collective_cause
Definition hdf5.h:962
std::string get_global_no_collective_cause()
Definition hdf5.cc:273
void write_none()
Definition hdf5.h:2190
Container read()
Definition hdf5.h:1766
bool get_query_io_mode() const
Definition hdf5.cc:296
std::uint32_t global_no_collective_cause
Definition hdf5.h:969
std::uint32_t get_global_no_collective_cause_as_hdf5_type()
Definition hdf5.cc:285
Container read_selection(const std::vector< hsize_t > &coordinates)
Definition hdf5.h:1801
H5D_mpio_actual_io_mode_t get_io_mode_as_hdf5_type()
Definition hdf5.cc:238
void set_query_io_mode(const bool new_query_io_mode)
Definition hdf5.cc:188
Container read_hyperslab(const std::vector< hsize_t > &offset, const std::vector< hsize_t > &count)
Definition hdf5.h:1853
unsigned int size
Definition hdf5.h:939
H5D_mpio_actual_io_mode_t io_mode
Definition hdf5.h:955
File(const std::string &name, const FileAccessMode mode)
Definition hdf5.cc:385
FileAccessMode
Definition hdf5.h:1079
DataSet open_dataset(const std::string &name) const
Definition hdf5.cc:378
Group open_group(const std::string &name) const
Definition hdf5.cc:362
Group(const std::string &name, const Group &parent_group, const bool mpi, const GroupAccessMode mode)
Definition hdf5.cc:319
DataSet create_dataset(const std::string &name, const std::vector< hsize_t > &dimensions) const
Definition hdf5.h:2230
Group create_group(const std::string &name) const
Definition hdf5.cc:370
GroupAccessMode
Definition hdf5.h:982
void write_dataset(const std::string &name, const Container &data) const
Definition hdf5.h:2241
std::shared_ptr< hid_t > hdf5_reference
Definition hdf5.h:413
std::string get_name() const
Definition hdf5.cc:76
T get_attribute(const std::string &attr_name) const
Definition hdf5.h:1594
const std::string name
Definition hdf5.h:404
void set_attribute(const std::string &attr_name, const T value)
Definition hdf5.h:1671
HDF5Object(const std::string &name, const bool mpi)
Definition hdf5.cc:68
const bool mpi
Definition hdf5.h:418
virtual size_type size() const override
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
const unsigned int DoFAccessor< structdim, dim, spacedim, level_dof_access >::dimension
#define DEAL_II_ASSERT_UNREACHABLE()
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)
const bool IsBlockVector< VectorType >::value
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:1494
std::string no_collective_cause_to_string(const std::uint32_t no_collective_cause)
Definition hdf5.h:1533
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:1418
unsigned int get_container_size(const std::vector< number > &data)
Definition hdf5.h:1389
void set_plist(hid_t &plist, const bool mpi)
Definition hdf5.h:1469
std::shared_ptr< hid_t > get_hdf5_datatype()
Definition hdf5.h:1274
std::vector< hsize_t > get_container_dimensions(const std::vector< number > &data)
Definition hdf5.h:1359
Definition hdf5.h:344