Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d5/d26/xdmf__utils_8h_source.html
DOLFINx 0.7.3
DOLFINx C++ interface
Loading...
Searching...
No Matches
xdmf_utils.h
1// Copyright (C) 2012 Chris N. Richardson
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{
29
30namespace fem
31{
32template <dolfinx::scalar T, std::floating_point U>
33class Function;
34} // namespace fem
35
36namespace fem
37{
38template <std::floating_point T>
39class CoordinateElement;
40class ElementDofLayout;
41} // namespace fem
42
43namespace mesh
44{
45template <std::floating_point T>
46class Mesh;
47class Topology;
48} // namespace mesh
49
50namespace io::xdmf_utils
51{
52
55std::pair<std::string, int> get_cell_type(const pugi::xml_node& topology_node);
56
59std::array<std::string, 2> get_hdf5_paths(const pugi::xml_node& dataitem_node);
60
61std::filesystem::path
62get_hdf5_filename(const std::filesystem::path& xdmf_filename);
63
65std::vector<std::int64_t> get_dataset_shape(const pugi::xml_node& dataset_node);
66
68std::int64_t get_num_cells(const pugi::xml_node& topology_node);
69
71std::string vtk_cell_type_str(mesh::CellType cell_type, int num_nodes);
72
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<
111 const std::int32_t,
112 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
113 xdofmap,
114 int entity_dim, std::span<const std::int64_t> entities,
115 std::span<const std::int32_t> data);
116
118template <typename T>
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)
123{
124 // Add DataItem node
125 assert(xml_node);
126 pugi::xml_node data_item_node = xml_node.append_child("DataItem");
127 assert(data_item_node);
128
129 // Add dimensions attribute
130 std::string dims;
131 for (auto d : shape)
132 dims += std::to_string(d) + std::string(" ");
133 dims.pop_back();
134 data_item_node.append_attribute("Dimensions") = dims.c_str();
135
136 // Set type for topology data (needed by XDMF to prevent default to
137 // float)
138 if (!number_type.empty())
139 data_item_node.append_attribute("NumberType") = number_type.c_str();
140
141 // Add format attribute
142 if (h5_id < 0)
143 {
144 data_item_node.append_attribute("Format") = "XML";
145 assert(shape.size() == 2);
146 std::ostringstream s;
147 s.precision(16);
148 for (std::size_t i = 0; i < x.size(); ++i)
149 {
150 if ((i + 1) % shape[1] == 0 and shape[1] != 0)
151 s << x.data()[i] << std::endl;
152 else
153 s << x.data()[i] << " ";
154 }
155
156 data_item_node.append_child(pugi::node_pcdata).set_value(s.str().c_str());
157 }
158 else
159 {
160 data_item_node.append_attribute("Format") = "HDF";
161
162 // Get name of HDF5 file, including path
163 const std::filesystem::path p = io::hdf5::get_filename(h5_id);
164 const std::filesystem::path filename = p.filename().c_str();
165
166 // Add HDF5 filename and HDF5 internal path to XML file
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());
170
171 // Compute data offset and range of values
172 std::int64_t local_shape0 = std::reduce(
173 std::next(shape.begin()), shape.end(), x.size(), std::divides{});
174
175 const std::array local_range{offset, offset + local_shape0};
176 io::hdf5::write_dataset(h5_id, h5_path, x.data(), local_range, shape,
177 use_mpi_io, false);
178
179 // Add partitioning attribute to dataset
180 // std::vector<std::size_t> partitions;
181 // std::vector<std::size_t> offset_tmp(1, offset);
182 // dolfinx::MPI::gather(comm, offset_tmp, partitions);
183 // dolfinx::MPI::broadcast(comm, partitions);
184 // io::hdf5::add_attribute(h5_id, h5_path, "partition", partitions);
185 }
186}
187
192template <typename T>
193std::vector<T> get_dataset(MPI_Comm comm, const pugi::xml_node& dataset_node,
194 hid_t h5_id,
195 std::array<std::int64_t, 2> range = {0, 0})
196{
197 // FIXME: Need to sort out datasset dimensions - can't depend on HDF5
198 // shape, and a Topology data item is not required to have a
199 // 'Dimensions' attribute since the dimensions can be determined from
200 // the number of cells and the cell type (for topology, one must
201 // supply cell type + (number of cells or dimensions).
202 //
203 // A geometry data item must have 'Dimensions' attribute.
204
205 assert(dataset_node);
206 pugi::xml_attribute format_attr = dataset_node.attribute("Format");
207 assert(format_attr);
208
209 // Get data set shape from 'Dimensions' attribute (empty if not
210 // available)
211 const std::vector shape_xml = xdmf_utils::get_dataset_shape(dataset_node);
212
213 const std::string format = format_attr.as_string();
214 std::vector<T> data_vector;
215 // Only read ASCII on process 0
216 const int mpi_rank = dolfinx::MPI::rank(comm);
217 if (format == "XML")
218 {
219 if (mpi_rank == 0)
220 {
221 // Read data and trim any leading/trailing whitespace
222 pugi::xml_node data_node = dataset_node.first_child();
223 assert(data_node);
224 std::string data_str = data_node.value();
225
226 // Split data based on spaces and line breaks
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"));
229
230 // Add data to numerical vector
231 data_vector.reserve(data_vector_str.size());
232 for (auto& v : data_vector_str)
233 {
234 if (v.begin() != v.end())
235 data_vector.push_back(
236 boost::lexical_cast<T>(boost::copy_range<std::string>(v)));
237 }
238 }
239 }
240 else if (format == "HDF")
241 {
242 // Get file and data path
243 auto paths = xdmf_utils::get_hdf5_paths(dataset_node);
244
245 // Get data shape from HDF5 file
246 const std::vector shape_hdf5 = io::hdf5::get_dataset_shape(h5_id, paths[1]);
247
248 // FIXME: should we support empty data sets?
249 // Check that data set is not empty
250 assert(!shape_hdf5.empty());
251 assert(shape_hdf5[0] != 0);
252
253 // Determine range of data to read from HDF5 file. This is
254 // complicated by the XML Dimension attribute and the HDF5 storage
255 // possibly having different shapes, e.g. the HDF5 storage may be a
256 // flat array.
257
258 // If range = {0, 0} then no range is supplied and we must determine
259 // the range
260 if (range[0] == 0 and range[1] == 0)
261 {
262 if (shape_xml == shape_hdf5)
263 {
264 range = dolfinx::MPI::local_range(mpi_rank, shape_hdf5[0],
265 dolfinx::MPI::size(comm));
266 }
267 else if (!shape_xml.empty() and shape_hdf5.size() == 1)
268 {
269 // Size of dims > 0
270 std::int64_t d = std::reduce(shape_xml.begin(), shape_xml.end(),
271 std::int64_t(1), std::multiplies{});
272
273 // Check for data size consistency
274 if (d * shape_xml[0] != shape_hdf5[0])
275 {
276 throw std::runtime_error("Data size in XDMF/XML and size of HDF5 "
277 "dataset are inconsistent");
278 }
279
280 // Compute data range to read
281 range = dolfinx::MPI::local_range(mpi_rank, shape_xml[0],
282 dolfinx::MPI::rank(comm));
283 range[0] *= d;
284 range[1] *= d;
285 }
286 else
287 {
288 throw std::runtime_error("This combination of array shapes in XDMF and "
289 "HDF5 is not supported");
290 }
291 }
292
293 // Retrieve data
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.");
297 else
298 {
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.");
302 }
303 }
304 else
305 throw std::runtime_error("Storage format \"" + format + "\" is unknown");
306
307 // Get dimensions for consistency (if available in DataItem node)
308 if (shape_xml.empty())
309 {
310 std::int64_t size = 1;
311 for (auto dim : shape_xml)
312 size *= dim;
313
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)
318 {
319 throw std::runtime_error(
320 "Data sizes in attribute and size of data read are inconsistent");
321 }
322 }
323
324 return data_vector;
325}
326
327} // namespace io::xdmf_utils
328} // namespace dolfinx
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