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);
106std::pair<std::vector<std::int32_t>, std::vector<std::int32_t>>
107distribute_entity_data(
108 const mesh::Topology& topology,
const std::vector<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>>
114 int entity_dim, std::span<const std::int64_t> entities,
115 std::span<const std::int32_t> data);
119void add_data_item(pugi::xml_node& xml_node, hid_t h5_id,
120 const std::string& h5_path, std::span<const T> x,
121 std::int64_t offset,
const std::vector<std::int64_t>& shape,
122 const std::string& number_type,
bool use_mpi_io)
126 pugi::xml_node data_item_node = xml_node.append_child(
"DataItem");
127 assert(data_item_node);
132 dims += std::
to_string(d) + std::string(
" ");
134 data_item_node.append_attribute(
"Dimensions") = dims.c_str();
138 if (!number_type.empty())
139 data_item_node.append_attribute(
"NumberType") = number_type.c_str();
144 data_item_node.append_attribute(
"Format") =
"XML";
145 assert(shape.size() == 2);
146 std::ostringstream s;
148 for (std::size_t i = 0; i < x.size(); ++i)
150 if ((i + 1) % shape[1] == 0 and shape[1] != 0)
151 s << x.data()[i] << std::endl;
153 s << x.data()[i] <<
" ";
156 data_item_node.append_child(pugi::node_pcdata).set_value(s.str().c_str());
160 data_item_node.append_attribute(
"Format") =
"HDF";
163 const std::filesystem::path p = io::hdf5::get_filename(h5_id);
164 const std::filesystem::path filename = p.filename().c_str();
167 const std::string xdmf_path
168 = filename.string() + std::string(
":") + h5_path;
169 data_item_node.append_child(pugi::node_pcdata).set_value(xdmf_path.c_str());
172 std::int64_t local_shape0 = std::reduce(
173 std::next(shape.begin()), shape.end(), x.size(), std::divides{});
175 const std::array
local_range{offset, offset + local_shape0};
176 io::hdf5::write_dataset(h5_id, h5_path, x.data(), local_range, shape,
193std::vector<T> get_dataset(MPI_Comm comm,
const pugi::xml_node& dataset_node,
195 std::array<std::int64_t, 2> range = {0, 0})
205 assert(dataset_node);
206 pugi::xml_attribute format_attr = dataset_node.attribute(
"Format");
211 const std::vector shape_xml = xdmf_utils::get_dataset_shape(dataset_node);
213 const std::string format = format_attr.as_string();
214 std::vector<T> data_vector;
222 pugi::xml_node data_node = dataset_node.first_child();
224 std::string data_str = data_node.value();
227 std::vector<boost::iterator_range<std::string::iterator>> data_vector_str;
228 boost::split(data_vector_str, data_str, boost::is_any_of(
" \n"));
231 data_vector.reserve(data_vector_str.size());
232 for (
auto& v : data_vector_str)
234 if (v.begin() != v.end())
235 data_vector.push_back(
236 boost::lexical_cast<T>(boost::copy_range<std::string>(v)));
240 else if (format ==
"HDF")
243 auto paths = xdmf_utils::get_hdf5_paths(dataset_node);
246 const std::vector shape_hdf5 = io::hdf5::get_dataset_shape(h5_id, paths[1]);
250 assert(!shape_hdf5.empty());
251 assert(shape_hdf5[0] != 0);
260 if (range[0] == 0 and range[1] == 0)
262 if (shape_xml == shape_hdf5)
267 else if (!shape_xml.empty() and shape_hdf5.size() == 1)
270 std::int64_t d = std::reduce(shape_xml.begin(), shape_xml.end(),
271 std::int64_t(1), std::multiplies{});
274 if (d * shape_xml[0] != shape_hdf5[0])
276 throw std::runtime_error(
"Data size in XDMF/XML and size of HDF5 "
277 "dataset are inconsistent");
288 throw std::runtime_error(
"This combination of array shapes in XDMF and "
289 "HDF5 is not supported");
294 if (hid_t dset_id = io::hdf5::open_dataset(h5_id, paths[1]);
295 dset_id == H5I_INVALID_HID)
296 throw std::runtime_error(
"Failed to open HDF5 global dataset.");
299 data_vector = io::hdf5::read_dataset<T>(dset_id, range,
true);
300 if (herr_t err = H5Dclose(dset_id); err < 0)
301 throw std::runtime_error(
"Failed to close HDF5 global dataset.");
305 throw std::runtime_error(
"Storage format \"" + format +
"\" is unknown");
308 if (shape_xml.empty())
310 std::int64_t
size = 1;
311 for (
auto dim : shape_xml)
314 std::int64_t size_global = 0;
315 const std::int64_t size_local = data_vector.size();
316 MPI_Allreduce(&size_local, &size_global, 1, MPI_INT64_T, MPI_SUM, comm);
317 if (size != size_global)
319 throw std::runtime_error(
320 "Data sizes in attribute and size of data read are inconsistent");
int size(MPI_Comm comm)
Return size of the group (number of processes) associated with the communicator.
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
std::string to_string(CellType type)
Get the cell string type for a cell type.
Definition cell_types.cpp:17
Top-level namespace.
Definition defines.h:12