File: XDMFFile.h

package info (click to toggle)
fenics-dolfinx 1%3A0.9.0-11
  • links: PTS, VCS
  • area: main
  • in suites: forky
  • size: 5,376 kB
  • sloc: cpp: 33,701; python: 22,338; makefile: 230; sh: 171; xml: 55
file content (212 lines) | stat: -rw-r--r-- 6,979 bytes parent folder | download | duplicates (4)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
// Copyright (C) 2012-2020 Chris N. Richardson, Garth N. Wells and Michal Habera
//
// This file is part of DOLFINx (https://www.fenicsproject.org)
//
// SPDX-License-Identifier:    LGPL-3.0-or-later

#pragma once

#include "HDF5Interface.h"
#include <concepts>
#include <dolfinx/common/MPI.h>
#include <dolfinx/mesh/cell_types.h>
#include <filesystem>
#include <memory>
#include <string>
#include <variant>

namespace pugi
{
class xml_node;
class xml_document;
} // namespace pugi

namespace dolfinx::fem
{
template <std::floating_point T>
class CoordinateElement;
template <dolfinx::scalar T, std::floating_point U>
class Function;
} // namespace dolfinx::fem

namespace dolfinx::mesh
{
template <std::floating_point T>
class Geometry;
enum class GhostMode : int;
template <std::floating_point T>
class Mesh;
template <typename T>
class MeshTags;
} // namespace dolfinx::mesh

namespace dolfinx::io
{

/// Read and write mesh::Mesh, fem::Function and other objects in
/// XDMF.
///
/// This class supports the output of meshes and functions in XDMF
/// (http://www.xdmf.org) format. It creates an XML file that describes
/// the data and points to a HDF5 file that stores the actual problem
/// data. Output of data in parallel is supported.
///
/// XDMF is not suitable for higher order geometries, as their currently
/// only supports 1st and 2nd order geometries.
class XDMFFile
{
public:
  /// File encoding type
  enum class Encoding
  {
    HDF5,
    ASCII
  };

  /// Default encoding type
  static const Encoding default_encoding = Encoding::HDF5;

  /// Constructor
  XDMFFile(MPI_Comm comm, const std::filesystem::path& filename,
           std::string file_mode, Encoding encoding = default_encoding);

  /// Move constructor
  XDMFFile(XDMFFile&&) = default;

  /// Destructor
  ~XDMFFile();

  /// Close the file
  ///
  /// This closes open underlying HDF5 file. In ASCII mode the XML file
  /// is closed each time it is written to or read from, so close() has
  /// no effect.
  void close();

  /// Save Mesh
  /// @param[in] mesh
  /// @param[in] xpath XPath where Mesh Grid will be written
  template <std::floating_point U>
  void write_mesh(const mesh::Mesh<U>& mesh,
                  std::string xpath = "/Xdmf/Domain");

  /// Save Geometry
  /// @param[in] geometry
  /// @param[in] name
  /// @param[in] xpath XPath of a node where Geometry will be inserted
  void write_geometry(const mesh::Geometry<double>& geometry, std::string name,
                      std::string xpath = "/Xdmf/Domain");

  /// Read in Mesh
  /// @param[in] element Element that describes the geometry of a cell
  /// @param[in] mode The type of ghosting/halo to use for the mesh when
  ///   distributed in parallel
  /// @param[in] name
  /// @param[in] xpath XPath where Mesh Grid is located
  /// @return A Mesh distributed on the same communicator as the
  ///   XDMFFile
  mesh::Mesh<double> read_mesh(const fem::CoordinateElement<double>& element,
                               mesh::GhostMode mode, std::string name,
                               std::string xpath = "/Xdmf/Domain") const;

  /// Read Topology data for Mesh
  /// @param[in] name Name of the mesh (Grid)
  /// @param[in] xpath XPath where Mesh Grid data is located
  /// @return (Cell type, degree), and cells topology (global node indexing)
  std::pair<std::vector<std::int64_t>, std::array<std::size_t, 2>>
  read_topology_data(std::string name,
                     std::string xpath = "/Xdmf/Domain") const;

  /// Read Geometry data for Mesh
  /// @param[in] name Name of the mesh (Grid)
  /// @param[in] xpath XPath where Mesh Grid data is located
  /// @return points on each process
  std::pair<std::variant<std::vector<float>, std::vector<double>>,
            std::array<std::size_t, 2>>
  read_geometry_data(std::string name,
                     std::string xpath = "/Xdmf/Domain") const;

  /// Read information about cell type
  /// @param[in] grid_name Name of Grid for which cell type is needed
  /// @param[in] xpath XPath where Grid is stored
  std::pair<mesh::CellType, int>
  read_cell_type(std::string grid_name, std::string xpath = "/Xdmf/Domain");

  /// @brief Write a fem::Function to file.
  ///
  /// @pre The fem::Function `u` must be (i) a lowest-order (P0)
  /// discontinuous Lagrange element or (ii) a continuous Lagrange
  /// element where the element 'nodes' are the same as the nodes of its
  /// mesh::Mesh. Otherwise an exception is raised.
  ///
  /// @note User interpolation to a suitable Lagrange space may be
  /// required to satisfy the precondition on `u`. The VTX output
  /// (io::VTXWriter) format is recommended over XDMF for discontinuous
  /// and/or high-order spaces.
  ///
  /// @param[in] u Function to write to file.
  /// @param[in] t Time stamp to associate with `u`.
  /// @param[in] mesh_xpath XPath for a Grid under which `u` will be
  /// inserted.
  template <dolfinx::scalar T, std::floating_point U = scalar_value_type_t<T>>
  void write_function(const fem::Function<T, U>& u, double t,
                      std::string mesh_xpath
                      = "/Xdmf/Domain/Grid[@GridType='Uniform'][1]");

  /// Write MeshTags
  /// @param[in] meshtags
  /// @param[in] x Mesh geometry
  /// @param[in] geometry_xpath XPath where Geometry is already stored
  /// in file
  /// @param[in] xpath XPath where MeshTags Grid will be inserted
  template <std::floating_point T>
  void write_meshtags(const mesh::MeshTags<std::int32_t>& meshtags,
                      const mesh::Geometry<T>& x, std::string geometry_xpath,
                      std::string xpath = "/Xdmf/Domain");

  /// Read MeshTags
  /// @param[in] mesh The Mesh that the data is defined on
  /// @param[in] name
  /// @param[in] xpath XPath where MeshTags Grid is stored in file
  mesh::MeshTags<std::int32_t>
  read_meshtags(const mesh::Mesh<double>& mesh, std::string name,
                std::string xpath = "/Xdmf/Domain");

  /// Write Information
  /// @param[in] name
  /// @param[in] value String to store into Information tag
  /// @param[in] xpath XPath where Information will be inserted
  void write_information(std::string name, std::string value,
                         std::string xpath = "/Xdmf/Domain/");

  /// Read Information
  /// @param[in] name
  /// @param[in] xpath XPath where Information is stored in file
  std::string read_information(std::string name,
                               std::string xpath = "/Xdmf/Domain/");

  /// Get the MPI communicator
  /// @return The MPI communicator for the file object
  MPI_Comm comm() const;

private:
  // MPI communicator
  dolfinx::MPI::Comm _comm;

  // Filename
  std::filesystem::path _filename;

  // File mode
  std::string _file_mode;

  // HDF5 file handle
  hid_t _h5_id;

  // The XML document currently representing the XDMF which needs to be
  // kept open for time series etc.
  std::unique_ptr<pugi::xml_document> _xml_doc;

  Encoding _encoding;
};

} // namespace dolfinx::io