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