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