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.5.1
DOLFINx C++ interface
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 <dolfinx/common/utils.h>
12 #include <dolfinx/mesh/cell_types.h>
13 #include <filesystem>
14 #include <numeric>
15 #include <pugixml.hpp>
16 #include <span>
17 #include <string>
18 #include <utility>
19 #include <vector>
20 
21 namespace pugi
22 {
23 class xml_node;
24 } // namespace pugi
25 
26 namespace dolfinx
27 {
28 
29 namespace fem
30 {
31 template <typename T>
32 class Function;
33 } // namespace fem
34 
35 namespace fem
36 {
37 class CoordinateElement;
38 }
39 
40 namespace mesh
41 {
42 class Mesh;
43 }
44 
45 namespace io::xdmf_utils
46 {
47 
48 // Get DOLFINx cell type string from XML topology node
49 // @return DOLFINx cell type and polynomial degree
50 std::pair<std::string, int> get_cell_type(const pugi::xml_node& topology_node);
51 
52 // Return (0) HDF5 filename and (1) path in HDF5 file from a DataItem
53 // node
54 std::array<std::string, 2> get_hdf5_paths(const pugi::xml_node& dataitem_node);
55 
56 std::filesystem::path
57 get_hdf5_filename(const std::filesystem::path& xdmf_filename);
58 
60 std::vector<std::int64_t> get_dataset_shape(const pugi::xml_node& dataset_node);
61 
63 std::int64_t get_num_cells(const pugi::xml_node& topology_node);
64 
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);
70 
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);
75 
77 std::string vtk_cell_type_str(mesh::CellType cell_type, int num_nodes);
78 
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);
108 
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)
115 {
116  // Add DataItem node
117  assert(xml_node);
118  pugi::xml_node data_item_node = xml_node.append_child("DataItem");
119  assert(data_item_node);
120 
121  // Add dimensions attribute
122  std::string dims;
123  for (auto d : shape)
124  dims += std::to_string(d) + std::string(" ");
125  dims.pop_back();
126  data_item_node.append_attribute("Dimensions") = dims.c_str();
127 
128  // Set type for topology data (needed by XDMF to prevent default to
129  // float)
130  if (!number_type.empty())
131  data_item_node.append_attribute("NumberType") = number_type.c_str();
132 
133  // Add format attribute
134  if (h5_id < 0)
135  {
136  data_item_node.append_attribute("Format") = "XML";
137  assert(shape.size() == 2);
138  std::ostringstream s;
139  s.precision(16);
140  for (std::size_t i = 0; i < (std::size_t)x.size(); ++i)
141  {
142  if ((i + 1) % shape[1] == 0 and shape[1] != 0)
143  s << x.data()[i] << std::endl;
144  else
145  s << x.data()[i] << " ";
146  }
147 
148  data_item_node.append_child(pugi::node_pcdata).set_value(s.str().c_str());
149  }
150  else
151  {
152  data_item_node.append_attribute("Format") = "HDF";
153 
154  // Get name of HDF5 file, including path
155  const std::filesystem::path p = HDF5Interface::get_filename(h5_id);
156  const std::filesystem::path filename = p.filename().c_str();
157 
158  // Add HDF5 filename and HDF5 internal path to XML file
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());
162 
163  // Compute data offset and range of values
164  std::int64_t local_shape0 = std::reduce(
165  std::next(shape.begin()), shape.end(), x.size(), std::divides{});
166 
167  const std::array local_range{offset, offset + local_shape0};
168  HDF5Interface::write_dataset(h5_id, h5_path, x.data(), local_range, shape,
169  use_mpi_io, false);
170 
171  // Add partitioning attribute to dataset
172  // std::vector<std::size_t> partitions;
173  // std::vector<std::size_t> offset_tmp(1, offset);
174  // dolfinx::MPI::gather(comm, offset_tmp, partitions);
175  // dolfinx::MPI::broadcast(comm, partitions);
176  // HDF5Interface::add_attribute(h5_id, h5_path, "partition", partitions);
177  }
178 }
179 
180 } // namespace io::xdmf_utils
181 } // namespace dolfinx
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