9 #include "HDF5Interface.h"
11 #include <dolfinx/common/utils.h>
12 #include <dolfinx/mesh/cell_types.h>
15 #include <pugixml.hpp>
37 class CoordinateElement;
45 namespace io::xdmf_utils
50 std::pair<std::string, int> get_cell_type(
const pugi::xml_node& topology_node);
54 std::array<std::string, 2> get_hdf5_paths(
const pugi::xml_node& dataitem_node);
57 get_hdf5_filename(
const std::filesystem::path& xdmf_filename);
60 std::vector<std::int64_t> get_dataset_shape(
const pugi::xml_node& dataset_node);
63 std::int64_t get_num_cells(
const pugi::xml_node& topology_node);
67 std::vector<double> get_point_data_values(
const fem::Function<double>& u);
68 std::vector<std::complex<double>>
69 get_point_data_values(
const fem::Function<std::complex<double>>& u);
72 std::vector<double> get_cell_data_values(
const fem::Function<double>& u);
73 std::vector<std::complex<double>>
74 get_cell_data_values(
const fem::Function<std::complex<double>>& u);
77 std::string vtk_cell_type_str(
mesh::CellType cell_type,
int num_nodes);
104 std::pair<std::vector<std::int32_t>, std::vector<std::int32_t>>
105 distribute_entity_data(
const mesh::Mesh& mesh,
int entity_dim,
106 const std::span<const std::int64_t>& entities,
107 const std::span<const std::int32_t>& data);
110 template <
typename T>
111 void add_data_item(pugi::xml_node& xml_node,
const hid_t h5_id,
112 const std::string& h5_path,
const T& x, std::int64_t offset,
113 const std::vector<std::int64_t>& shape,
114 const std::string& number_type,
bool use_mpi_io)
118 pugi::xml_node data_item_node = xml_node.append_child(
"DataItem");
119 assert(data_item_node);
124 dims += std::to_string(d) + std::string(
" ");
126 data_item_node.append_attribute(
"Dimensions") = dims.c_str();
130 if (!number_type.empty())
131 data_item_node.append_attribute(
"NumberType") = number_type.c_str();
136 data_item_node.append_attribute(
"Format") =
"XML";
137 assert(shape.size() == 2);
138 std::ostringstream s;
140 for (std::size_t i = 0; i < (std::size_t)x.size(); ++i)
142 if ((i + 1) % shape[1] == 0 and shape[1] != 0)
143 s << x.data()[i] << std::endl;
145 s << x.data()[i] <<
" ";
148 data_item_node.append_child(pugi::node_pcdata).set_value(s.str().c_str());
152 data_item_node.append_attribute(
"Format") =
"HDF";
156 const std::filesystem::path filename = p.filename().c_str();
159 const std::string xdmf_path
160 = filename.string() + std::string(
":") + h5_path;
161 data_item_node.append_child(pugi::node_pcdata).set_value(xdmf_path.c_str());
164 std::int64_t local_shape0 = std::reduce(
165 std::next(shape.begin()), shape.end(), x.size(), std::divides{});
167 const std::array
local_range{offset, offset + local_shape0};
static std::filesystem::path get_filename(hid_t handle)
Get filename.
Definition: HDF5Interface.cpp:136
static void write_dataset(const hid_t handle, const std::string &dataset_path, const T *data, const std::array< std::int64_t, 2 > &range, const std::vector< std::int64_t > &global_size, bool use_mpi_io, bool use_chunking)
Write data to existing HDF file as defined by range blocks on each process.
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:83
CellType
Cell type identifier.
Definition: cell_types.h:22