Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/v0.1.0/v0.9.0/cpp
DOLFINx  0.1.0
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 "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 #include <xtl/xspan.hpp>
19 
20 namespace pugi
21 {
22 class xml_node;
23 } // namespace pugi
24 
25 namespace dolfinx
26 {
27 
28 namespace fem
29 {
30 template <typename T>
31 class Function;
32 } // namespace fem
33 
34 namespace fem
35 {
36 class CoordinateElement;
37 }
38 
39 namespace mesh
40 {
41 class Mesh;
42 }
43 
44 namespace io::xdmf_utils
45 {
46 
47 // Get DOLFINx cell type string from XML topology node
48 // @return DOLFINx cell type and polynomial degree
49 std::pair<std::string, int> get_cell_type(const pugi::xml_node& topology_node);
50 
51 // Return (0) HDF5 filename and (1) path in HDF5 file from a DataItem
52 // node
53 std::array<std::string, 2> get_hdf5_paths(const pugi::xml_node& dataitem_node);
54 
55 std::string get_hdf5_filename(std::string xdmf_filename);
56 
58 std::vector<std::int64_t> get_dataset_shape(const pugi::xml_node& dataset_node);
59 
61 std::int64_t get_num_cells(const pugi::xml_node& topology_node);
62 
65 std::vector<double> get_point_data_values(const fem::Function<double>& u);
66 std::vector<std::complex<double>>
67 get_point_data_values(const fem::Function<std::complex<double>>& u);
68 
70 std::vector<double> get_cell_data_values(const fem::Function<double>& u);
71 std::vector<std::complex<double>>
72 get_cell_data_values(const fem::Function<std::complex<double>>& u);
73 
75 std::string vtk_cell_type_str(mesh::CellType cell_type, int num_nodes);
76 
96 std::pair<xt::xtensor<std::int32_t, 2>, std::vector<std::int32_t>>
97 extract_local_entities(const mesh::Mesh& mesh, int entity_dim,
98  const xt::xtensor<std::int64_t, 2>& entities,
99  const xtl::span<const std::int32_t>& values);
100 
102 template <typename T>
103 void add_data_item(pugi::xml_node& xml_node, const hid_t h5_id,
104  const std::string h5_path, const T& x,
105  const std::int64_t offset,
106  const std::vector<std::int64_t> shape,
107  const std::string number_type, const bool use_mpi_io)
108 {
109  // Add DataItem node
110  assert(xml_node);
111  pugi::xml_node data_item_node = xml_node.append_child("DataItem");
112  assert(data_item_node);
113 
114  // Add dimensions attribute
115  std::string dims;
116  for (auto d : shape)
117  dims += std::to_string(d) + " ";
118  dims.pop_back();
119  data_item_node.append_attribute("Dimensions") = dims.c_str();
120 
121  // Set type for topology data (needed by XDMF to prevent default to
122  // float)
123  if (!number_type.empty())
124  data_item_node.append_attribute("NumberType") = number_type.c_str();
125 
126  // Add format attribute
127  if (h5_id < 0)
128  {
129  data_item_node.append_attribute("Format") = "XML";
130  assert(shape.size() == 2);
131  data_item_node.append_child(pugi::node_pcdata)
132  .set_value(common::container_to_string(x, 16, shape[1]).c_str());
133  }
134  else
135  {
136  data_item_node.append_attribute("Format") = "HDF";
137 
138  // Get name of HDF5 file, including path
139  const std::string hdf5_filename = HDF5Interface::get_filename(h5_id);
140  const std::string filename = dolfinx::io::get_filename(hdf5_filename);
141 
142  // Add HDF5 filename and HDF5 internal path to XML file
143  const std::string xdmf_path = filename + ":" + h5_path;
144  data_item_node.append_child(pugi::node_pcdata).set_value(xdmf_path.c_str());
145 
146  // Compute total number of items and check for consistency with shape
147  assert(!shape.empty());
148  std::int64_t num_items_total = 1;
149  for (auto n : shape)
150  num_items_total *= n;
151 
152  // Compute data offset and range of values
153  std::int64_t local_shape0 = x.size();
154  for (std::size_t i = 1; i < shape.size(); ++i)
155  {
156  assert(local_shape0 % shape[i] == 0);
157  local_shape0 /= shape[i];
158  }
159 
160  const std::array local_range{offset, offset + local_shape0};
161  HDF5Interface::write_dataset(h5_id, h5_path, x.data(), local_range, shape,
162  use_mpi_io, false);
163 
164  // Add partitioning attribute to dataset
165  // std::vector<std::size_t> partitions;
166  // std::vector<std::size_t> offset_tmp(1, offset);
167  // dolfinx::MPI::gather(comm, offset_tmp, partitions);
168  // dolfinx::MPI::broadcast(comm, partitions);
169  // HDF5Interface::add_attribute(h5_id, h5_path, "partition", partitions);
170  }
171 }
172 
173 } // namespace io::xdmf_utils
174 } // namespace dolfinx
dolfinx::mesh::CellType
CellType
Cell type identifier.
Definition: cell_types.h:21
dolfinx::common::container_to_string
std::string container_to_string(const T &x, const int precision, const int linebreak)
Convert a container to string.
Definition: utils.h:62
dolfinx::io::get_filename
std::string get_filename(const std::string &fullname)
Get filename from a fully qualified path and filename.
Definition: utils.cpp:13
dolfinx::io::HDF5Interface::write_dataset
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.
dolfinx::io::HDF5Interface::get_filename
static std::string get_filename(hid_t handle)
Get filename.
Definition: HDF5Interface.cpp:128