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>
32template <dolfinx::scalar T, std::
floating_po
int U>
38template <std::
floating_po
int T>
39class CoordinateElement;
40class ElementDofLayout;
45template <std::
floating_po
int T>
50namespace io::xdmf_utils
55std::pair<std::string, int> get_cell_type(
const pugi::xml_node& topology_node);
59std::array<std::string, 2> get_hdf5_paths(
const pugi::xml_node& dataitem_node);
62get_hdf5_filename(
const std::filesystem::path& xdmf_filename);
65std::vector<std::int64_t> get_dataset_shape(
const pugi::xml_node& dataset_node);
68std::int64_t get_num_cells(
const pugi::xml_node& topology_node);
71std::string vtk_cell_type_str(
mesh::CellType cell_type,
int num_nodes);
108std::pair<std::vector<std::int32_t>, std::vector<T>> distribute_entity_data(
109 const mesh::Topology& topology, std::span<const std::int64_t> nodes_g,
110 std::int64_t num_nodes_g,
const fem::ElementDofLayout& cmap_dof_layout,
111 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
113 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
116 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
118 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
120 std::span<const T> data);
124void add_data_item(pugi::xml_node& xml_node, hid_t h5_id,
125 const std::string& h5_path, std::span<const T> x,
126 std::int64_t offset,
const std::vector<std::int64_t>& shape,
127 const std::string& number_type,
bool use_mpi_io)
131 pugi::xml_node data_item_node = xml_node.append_child(
"DataItem");
132 assert(data_item_node);
137 dims += std::to_string(d) + std::string(
" ");
139 data_item_node.append_attribute(
"Dimensions") = dims.c_str();
143 if (!number_type.empty())
144 data_item_node.append_attribute(
"NumberType") = number_type.c_str();
149 data_item_node.append_attribute(
"Format") =
"XML";
150 assert(shape.size() == 2);
151 std::ostringstream s;
153 for (std::size_t i = 0; i < x.size(); ++i)
155 if ((i + 1) % shape[1] == 0 and shape[1] != 0)
156 s << x.data()[i] << std::endl;
158 s << x.data()[i] <<
" ";
161 data_item_node.append_child(pugi::node_pcdata).set_value(s.str().c_str());
165 data_item_node.append_attribute(
"Format") =
"HDF";
168 const std::filesystem::path p = io::hdf5::get_filename(h5_id);
169 const std::filesystem::path filename = p.filename().c_str();
172 const std::string xdmf_path
173 = filename.string() + std::string(
":") + h5_path;
174 data_item_node.append_child(pugi::node_pcdata).set_value(xdmf_path.c_str());
177 std::int64_t local_shape0 = std::reduce(
178 std::next(shape.begin()), shape.end(), x.size(), std::divides{});
180 const std::array
local_range{offset, offset + local_shape0};
181 io::hdf5::write_dataset(h5_id, h5_path, x.data(), local_range, shape,
198std::vector<T> get_dataset(MPI_Comm comm,
const pugi::xml_node& dataset_node,
200 std::array<std::int64_t, 2> range = {0, 0})
210 assert(dataset_node);
211 pugi::xml_attribute format_attr = dataset_node.attribute(
"Format");
216 const std::vector shape_xml = xdmf_utils::get_dataset_shape(dataset_node);
218 const std::string format = format_attr.as_string();
219 std::vector<T> data_vector;
227 pugi::xml_node data_node = dataset_node.first_child();
229 std::string data_str = data_node.value();
232 std::vector<boost::iterator_range<std::string::iterator>> data_vector_str;
233 boost::split(data_vector_str, data_str, boost::is_any_of(
" \n"));
236 data_vector.reserve(data_vector_str.size());
237 for (
auto& v : data_vector_str)
239 if (v.begin() != v.end())
240 data_vector.push_back(
241 boost::lexical_cast<T>(boost::copy_range<std::string>(v)));
245 else if (format ==
"HDF")
248 auto paths = xdmf_utils::get_hdf5_paths(dataset_node);
251 const std::vector shape_hdf5 = io::hdf5::get_dataset_shape(h5_id, paths[1]);
255 assert(!shape_hdf5.empty());
256 assert(shape_hdf5[0] != 0);
265 if (range[0] == 0 and range[1] == 0)
267 if (shape_xml == shape_hdf5)
272 else if (!shape_xml.empty() and shape_hdf5.size() == 1)
275 std::int64_t d = std::reduce(shape_xml.begin(), shape_xml.end(),
276 std::int64_t(1), std::multiplies{});
279 if (d * shape_xml[0] != shape_hdf5[0])
281 throw std::runtime_error(
"Data size in XDMF/XML and size of HDF5 "
282 "dataset are inconsistent");
293 throw std::runtime_error(
"This combination of array shapes in XDMF and "
294 "HDF5 is not supported");
299 if (hid_t dset_id = io::hdf5::open_dataset(h5_id, paths[1]);
300 dset_id == H5I_INVALID_HID)
301 throw std::runtime_error(
"Failed to open HDF5 global dataset.");
304 data_vector = io::hdf5::read_dataset<T>(dset_id, range,
true);
305 if (herr_t err = H5Dclose(dset_id); err < 0)
306 throw std::runtime_error(
"Failed to close HDF5 global dataset.");
310 throw std::runtime_error(
"Storage format \"" + format +
"\" is unknown");
313 if (shape_xml.empty())
315 std::int64_t
size = 1;
316 for (
auto dim : shape_xml)
319 std::int64_t size_global = 0;
320 const std::int64_t size_local = data_vector.size();
321 MPI_Allreduce(&size_local, &size_global, 1, MPI_INT64_T, MPI_SUM, comm);
322 if (size != size_global)
324 throw std::runtime_error(
325 "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