9 #include "HDF5Interface.h"
10 #include "pugixml.hpp"
13 #include <dolfinx/common/utils.h>
14 #include <dolfinx/mesh/cell_types.h>
35 class CoordinateElement;
43 namespace io::xdmf_utils
48 std::pair<std::string, int> get_cell_type(
const pugi::xml_node& topology_node);
52 std::array<std::string, 2> get_hdf5_paths(
const pugi::xml_node& dataitem_node);
54 std::string get_hdf5_filename(std::string xdmf_filename);
57 std::vector<std::int64_t> get_dataset_shape(
const pugi::xml_node& dataset_node);
60 std::int64_t get_num_cells(
const pugi::xml_node& topology_node);
64 std::vector<double> get_point_data_values(
const fem::Function<double>& u);
65 std::vector<std::complex<double>>
66 get_point_data_values(
const fem::Function<std::complex<double>>& u);
69 std::vector<double> get_cell_data_values(
const fem::Function<double>& u);
70 std::vector<std::complex<double>>
71 get_cell_data_values(
const fem::Function<std::complex<double>>& u);
74 std::string vtk_cell_type_str(
mesh::CellType cell_type,
int num_nodes);
96 Eigen::Array<std::int32_t, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>,
97 std::vector<std::int32_t>>
98 extract_local_entities(
99 const mesh::Mesh& mesh,
const int entity_dim,
100 const Eigen::Array<std::int64_t, Eigen::Dynamic, Eigen::Dynamic,
101 Eigen::RowMajor>& entities,
102 const std::vector<std::int32_t>& values);
105 template <
typename T>
106 void add_data_item(pugi::xml_node& xml_node,
const hid_t h5_id,
107 const std::string h5_path,
const T& x,
108 const std::int64_t offset,
109 const std::vector<std::int64_t> shape,
110 const std::string number_type,
const bool use_mpi_io)
114 pugi::xml_node data_item_node = xml_node.append_child(
"DataItem");
115 assert(data_item_node);
120 dims += std::to_string(d) +
" ";
122 data_item_node.append_attribute(
"Dimensions") = dims.c_str();
126 if (!number_type.empty())
127 data_item_node.append_attribute(
"NumberType") = number_type.c_str();
132 data_item_node.append_attribute(
"Format") =
"XML";
133 assert(shape.size() == 2);
134 data_item_node.append_child(pugi::node_pcdata)
139 data_item_node.append_attribute(
"Format") =
"HDF";
146 const std::string xdmf_path = filename +
":" + h5_path;
147 data_item_node.append_child(pugi::node_pcdata).set_value(xdmf_path.c_str());
150 assert(!shape.empty());
151 std::int64_t num_items_total = 1;
153 num_items_total *= n;
156 std::int64_t local_shape0 = x.size();
157 for (std::size_t i = 1; i < shape.size(); ++i)
159 assert(local_shape0 % shape[i] == 0);
160 local_shape0 /= shape[i];
163 const std::array local_range{offset, offset + local_shape0};
static std::string get_filename(hid_t handle)
Get filename.
Definition: HDF5Interface.cpp:129
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.
std::string container_to_string(const T &x, const int precision, const int linebreak)
Convert a container to string.
Definition: utils.h:63
std::string get_filename(const std::string &fullname)
Get filename from a fully qualified path and filename.
Definition: utils.cpp:13
CellType
Cell type identifier.
Definition: cell_types.h:21