DOLFINx 0.10.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
xdmf_utils.h
1// Copyright (C) 2012-2024 Chris N. Richardson, Jørgen S. Dokken
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#include "HDF5Interface.h"
10#include <array>
11#include <basix/mdspan.hpp>
12#include <boost/algorithm/string.hpp>
13#include <boost/lexical_cast.hpp>
14#include <complex>
15#include <concepts>
16#include <dolfinx/common/MPI.h>
17#include <dolfinx/common/types.h>
18#include <dolfinx/mesh/cell_types.h>
19#include <filesystem>
20#include <numeric>
21#include <pugixml.hpp>
22#include <span>
23#include <string>
24#include <utility>
25#include <vector>
26
27namespace dolfinx
28{
29namespace fem
30{
31template <dolfinx::scalar T, std::floating_point U>
32class Function;
33} // namespace fem
34
35namespace fem
36{
37template <std::floating_point T>
38class CoordinateElement;
39class ElementDofLayout;
40} // namespace fem
41
42namespace mesh
43{
44template <std::floating_point T>
45class Mesh;
46class Topology;
47} // namespace mesh
48
49namespace io::xdmf_utils
50{
51
54std::pair<std::string, int> get_cell_type(const pugi::xml_node& topology_node);
55
58std::array<std::string, 2> get_hdf5_paths(const pugi::xml_node& dataitem_node);
59
60std::filesystem::path
61get_hdf5_filename(const std::filesystem::path& xdmf_filename);
62
64std::vector<std::int64_t> get_dataset_shape(const pugi::xml_node& dataset_node);
65
67std::int64_t get_num_cells(const pugi::xml_node& topology_node);
68
70std::string vtk_cell_type_str(mesh::CellType cell_type, int num_nodes);
71
106template <typename T>
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<
111 const std::int32_t,
112 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
113 xdofmap,
114 int entity_dim,
115 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
116 const std::int64_t,
117 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
118 entities,
119 std::span<const T> data);
120
122template <typename T>
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)
127{
128 // Add DataItem node
129 assert(xml_node);
130 pugi::xml_node data_item_node = xml_node.append_child("DataItem");
131 assert(data_item_node);
132
133 // Add dimensions attribute
134 std::string dims;
135 for (auto d : shape)
136 dims += std::to_string(d) + std::string(" ");
137 dims.pop_back();
138 data_item_node.append_attribute("Dimensions") = dims.c_str();
139
140 // Set type for topology data (needed by XDMF to prevent default to
141 // float)
142 if (!number_type.empty())
143 data_item_node.append_attribute("NumberType") = number_type.c_str();
144
145 // Add format attribute
146 if (h5_id < 0)
147 {
148 data_item_node.append_attribute("Format") = "XML";
149 assert(shape.size() == 2);
150 std::ostringstream s;
151 s.precision(16);
152 for (std::size_t i = 0; i < x.size(); ++i)
153 {
154 if ((i + 1) % shape[1] == 0 and shape[1] != 0)
155 s << x.data()[i] << std::endl;
156 else
157 s << x.data()[i] << " ";
158 }
159
160 data_item_node.append_child(pugi::node_pcdata).set_value(s.str().c_str());
161 }
162 else
163 {
164 data_item_node.append_attribute("Format") = "HDF";
165
166 // Get name of HDF5 file, including path
167 const std::filesystem::path p = io::hdf5::get_filename(h5_id);
168 const std::filesystem::path filename = p.filename().c_str();
169
170 // Add HDF5 filename and HDF5 internal path to XML file
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());
174
175 // Compute data offset and range of values
176 std::int64_t local_shape0 = std::reduce(
177 std::next(shape.begin()), shape.end(), x.size(), std::divides{});
178
179 const std::array local_range{offset, offset + local_shape0};
180 io::hdf5::write_dataset(h5_id, h5_path, x.data(), local_range, shape,
181 use_mpi_io, false);
182
183 // Add partitioning attribute to dataset
184 // std::vector<std::size_t> partitions;
185 // std::vector<std::size_t> offset_tmp(1, offset);
186 // dolfinx::MPI::gather(comm, offset_tmp, partitions);
187 // dolfinx::MPI::broadcast(comm, partitions);
188 // io::hdf5::add_attribute(h5_id, h5_path, "partition", partitions);
189 }
190}
191
196template <typename T>
197std::vector<T> get_dataset(MPI_Comm comm, const pugi::xml_node& dataset_node,
198 hid_t h5_id,
199 std::array<std::int64_t, 2> range = {0, 0})
200{
201 // FIXME: Need to sort out dataset dimensions - can't depend on HDF5
202 // shape, and a Topology data item is not required to have a
203 // 'Dimensions' attribute since the dimensions can be determined from
204 // the number of cells and the cell type (for topology, one must
205 // supply cell type + (number of cells or dimensions)).
206 //
207 // A geometry data item must have 'Dimensions' attribute.
208
209 assert(dataset_node);
210 pugi::xml_attribute format_attr = dataset_node.attribute("Format");
211 assert(format_attr);
212
213 // Get data set shape from 'Dimensions' attribute (empty if not
214 // available)
215 const std::vector shape_xml = xdmf_utils::get_dataset_shape(dataset_node);
216
217 const std::string format = format_attr.as_string();
218 std::vector<T> data_vector;
219 // Only read ASCII on process 0
220 const int mpi_rank = dolfinx::MPI::rank(comm);
221 if (format == "XML")
222 {
223 if (mpi_rank == 0)
224 {
225 // Read data and trim any leading/trailing whitespace
226 pugi::xml_node data_node = dataset_node.first_child();
227 assert(data_node);
228 std::string data_str = data_node.value();
229
230 // Split data based on spaces and line breaks
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"));
233
234 // Add data to numerical vector
235 data_vector.reserve(data_vector_str.size());
236 for (auto& v : data_vector_str)
237 {
238 if (v.begin() != v.end())
239 data_vector.push_back(
240 boost::lexical_cast<T>(boost::copy_range<std::string>(v)));
241 }
242 }
243 }
244 else if (format == "HDF")
245 {
246 // Get file and data path
247 auto paths = xdmf_utils::get_hdf5_paths(dataset_node);
248
249 // Get data shape from HDF5 file
250 const std::vector shape_hdf5 = io::hdf5::get_dataset_shape(h5_id, paths[1]);
251
252 // FIXME: should we support empty data sets?
253 // Check that data set is not empty
254 assert(!shape_hdf5.empty());
255 assert(shape_hdf5[0] != 0);
256
257 // Determine range of data to read from HDF5 file. This is
258 // complicated by the XML Dimension attribute and the HDF5 storage
259 // possibly having different shapes, e.g. the HDF5 storage may be a
260 // flat array.
261
262 // If range = {0, 0} then no range is supplied and we must determine
263 // the range
264 if (range[0] == 0 and range[1] == 0)
265 {
266 if (shape_xml == shape_hdf5)
267 {
268 range = dolfinx::MPI::local_range(mpi_rank, shape_hdf5[0],
269 dolfinx::MPI::size(comm));
270 }
271 else if (!shape_xml.empty() and shape_hdf5.size() == 1)
272 {
273 // Size of dims > 0
274 std::int64_t d = std::reduce(shape_xml.begin(), shape_xml.end(),
275 std::int64_t(1), std::multiplies{});
276
277 // Check for data size consistency
278 if (d * shape_xml[0] != shape_hdf5[0])
279 {
280 throw std::runtime_error("Data size in XDMF/XML and size of HDF5 "
281 "dataset are inconsistent");
282 }
283
284 // Compute data range to read
285 range = dolfinx::MPI::local_range(mpi_rank, shape_xml[0],
286 dolfinx::MPI::rank(comm));
287 range[0] *= d;
288 range[1] *= d;
289 }
290 else
291 {
292 throw std::runtime_error("This combination of array shapes in XDMF and "
293 "HDF5 is not supported");
294 }
295 }
296
297 // Retrieve data
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.");
301 else
302 {
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.");
306 }
307 }
308 else
309 throw std::runtime_error("Storage format \"" + format + "\" is unknown");
310
311 // Get dimensions for consistency (if available in DataItem node)
312 if (shape_xml.empty())
313 {
314 std::int64_t size = 1;
315 for (auto dim : shape_xml)
316 size *= dim;
317
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)
322 {
323 throw std::runtime_error(
324 "Data sizes in attribute and size of data read are inconsistent");
325 }
326 }
327
328 return data_vector;
329}
330
331} // namespace io::xdmf_utils
332} // namespace dolfinx
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