9#include "HDF5Interface.h"
11#include <basix/mdspan.hpp>
12#include <boost/algorithm/string.hpp>
13#include <boost/lexical_cast.hpp>
16#include <dolfinx/common/MPI.h>
17#include <dolfinx/common/types.h>
18#include <dolfinx/mesh/cell_types.h>
31template <dolfinx::scalar T, std::
floating_po
int U>
37template <std::
floating_po
int T>
38class CoordinateElement;
39class ElementDofLayout;
44template <std::
floating_po
int T>
49namespace io::xdmf_utils
54std::pair<std::string, int> get_cell_type(
const pugi::xml_node& topology_node);
58std::array<std::string, 2> get_hdf5_paths(
const pugi::xml_node& dataitem_node);
61get_hdf5_filename(
const std::filesystem::path& xdmf_filename);
64std::vector<std::int64_t> get_dataset_shape(
const pugi::xml_node& dataset_node);
67std::int64_t get_num_cells(
const pugi::xml_node& topology_node);
70std::string vtk_cell_type_str(
mesh::CellType cell_type,
int num_nodes);
107std::pair<std::vector<std::int32_t>, std::vector<T>> distribute_entity_data(
108 const mesh::Topology& topology, std::span<const std::int64_t> nodes_g,
109 std::int64_t num_nodes_g,
const fem::ElementDofLayout& cmap_dof_layout,
110 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
112 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
115 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
117 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
119 std::span<const T> data);
123void add_data_item(pugi::xml_node& xml_node, hid_t h5_id,
124 const std::string& h5_path, std::span<const T> x,
125 std::int64_t offset,
const std::vector<std::int64_t>& shape,
126 const std::string& number_type,
bool use_mpi_io)
130 pugi::xml_node data_item_node = xml_node.append_child(
"DataItem");
131 assert(data_item_node);
136 dims += std::to_string(d) + std::string(
" ");
138 data_item_node.append_attribute(
"Dimensions") = dims.c_str();
142 if (!number_type.empty())
143 data_item_node.append_attribute(
"NumberType") = number_type.c_str();
148 data_item_node.append_attribute(
"Format") =
"XML";
149 assert(shape.size() == 2);
150 std::ostringstream s;
152 for (std::size_t i = 0; i < x.size(); ++i)
154 if ((i + 1) % shape[1] == 0 and shape[1] != 0)
155 s << x.data()[i] << std::endl;
157 s << x.data()[i] <<
" ";
160 data_item_node.append_child(pugi::node_pcdata).set_value(s.str().c_str());
164 data_item_node.append_attribute(
"Format") =
"HDF";
167 const std::filesystem::path p = io::hdf5::get_filename(h5_id);
168 const std::filesystem::path filename = p.filename().c_str();
171 const std::string xdmf_path
172 = filename.string() + std::string(
":") + h5_path;
173 data_item_node.append_child(pugi::node_pcdata).set_value(xdmf_path.c_str());
176 std::int64_t local_shape0 = std::reduce(
177 std::next(shape.begin()), shape.end(), x.size(), std::divides{});
179 const std::array
local_range{offset, offset + local_shape0};
180 io::hdf5::write_dataset(h5_id, h5_path, x.data(), local_range, shape,
197std::vector<T> get_dataset(MPI_Comm comm,
const pugi::xml_node& dataset_node,
199 std::array<std::int64_t, 2> range = {0, 0})
209 assert(dataset_node);
210 pugi::xml_attribute format_attr = dataset_node.attribute(
"Format");
215 const std::vector shape_xml = xdmf_utils::get_dataset_shape(dataset_node);
217 const std::string format = format_attr.as_string();
218 std::vector<T> data_vector;
226 pugi::xml_node data_node = dataset_node.first_child();
228 std::string data_str = data_node.value();
231 std::vector<boost::iterator_range<std::string::iterator>> data_vector_str;
232 boost::split(data_vector_str, data_str, boost::is_any_of(
" \n"));
235 data_vector.reserve(data_vector_str.size());
236 for (
auto& v : data_vector_str)
238 if (v.begin() != v.end())
239 data_vector.push_back(
240 boost::lexical_cast<T>(boost::copy_range<std::string>(v)));
244 else if (format ==
"HDF")
247 auto paths = xdmf_utils::get_hdf5_paths(dataset_node);
250 const std::vector shape_hdf5 = io::hdf5::get_dataset_shape(h5_id, paths[1]);
254 assert(!shape_hdf5.empty());
255 assert(shape_hdf5[0] != 0);
264 if (range[0] == 0 and range[1] == 0)
266 if (shape_xml == shape_hdf5)
271 else if (!shape_xml.empty() and shape_hdf5.size() == 1)
274 std::int64_t d = std::reduce(shape_xml.begin(), shape_xml.end(),
275 std::int64_t(1), std::multiplies{});
278 if (d * shape_xml[0] != shape_hdf5[0])
280 throw std::runtime_error(
"Data size in XDMF/XML and size of HDF5 "
281 "dataset are inconsistent");
292 throw std::runtime_error(
"This combination of array shapes in XDMF and "
293 "HDF5 is not supported");
298 if (hid_t dset_id = io::hdf5::open_dataset(h5_id, paths[1]);
299 dset_id == H5I_INVALID_HID)
300 throw std::runtime_error(
"Failed to open HDF5 global dataset.");
303 data_vector = io::hdf5::read_dataset<T>(dset_id, range,
true);
304 if (herr_t err = H5Dclose(dset_id); err < 0)
305 throw std::runtime_error(
"Failed to close HDF5 global dataset.");
309 throw std::runtime_error(
"Storage format \"" + format +
"\" is unknown");
312 if (shape_xml.empty())
314 std::int64_t
size = 1;
315 for (
auto dim : shape_xml)
318 std::int64_t size_global = 0;
319 const std::int64_t size_local = data_vector.size();
320 MPI_Allreduce(&size_local, &size_global, 1, MPI_INT64_T, MPI_SUM, comm);
321 if (size != size_global)
323 throw std::runtime_error(
324 "Data sizes in attribute and size of data read are inconsistent");
int size(MPI_Comm comm)
Definition MPI.cpp:72
int rank(MPI_Comm comm)
Return process rank for the communicator.
Definition MPI.cpp:64
constexpr std::array< std::int64_t, 2 > local_range(int rank, std::int64_t N, int size)
Return local range for the calling process, partitioning the global [0, N - 1] range across all ranks...
Definition MPI.h:91
CellType
Cell type identifier.
Definition cell_types.h:22
Top-level namespace.
Definition defines.h:12