MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
vtkhdf.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#ifndef MFEM_VTKHDF
13#define MFEM_VTKHDF
14
15#include "../config/config.hpp"
16
17#ifdef MFEM_USE_HDF5
18
19#include "../fem/gridfunc.hpp"
20
21#include <hdf5.h>
22#include <cstdint>
23#include <unordered_map>
24
25#if defined(MFEM_USE_MPI) && defined(H5_HAVE_PARALLEL)
26#define MFEM_PARALLEL_HDF5
27#endif
28
29namespace mfem
30{
31
32/// @brief Low-level class for writing %VTKHDF data (for use in ParaView).
33///
34/// Users should typically use ParaViewHDFDataCollection, Mesh::SaveVTKHDF, or
35/// GridFunction::SaveVTKHDF instead.
36class VTKHDF
37{
38public:
39 /// Helper used in VTKHDF::VTKHDF for enabling/disabling restart mode.
40 struct Restart
41 {
42 bool enabled = false;
43 real_t time = 0.0;
44 Restart() = default;
45 Restart(real_t time_) : enabled(true), time(time_) { }
46 Restart(bool enabled_, real_t time_) : enabled(enabled_), time(time_) { }
47 static Restart Disabled() { return Restart(); }
48 };
49
50private:
51#ifdef MFEM_USE_MPI
52 /// MPI communicator (only available if MPI is enabled).
53 MPI_Comm comm = MPI_COMM_NULL;
54#endif
55 /// Size of the MPI communicator (1 if MPI is not enabled).
56 const int mpi_size = 1;
57 /// Rank within MPI communicator (0 if MPI is not enabled).
58 const int mpi_rank = 0;
59
60 /// File access property list (needed for MPI I/O).
61 hid_t fapl = H5I_INVALID_HID;
62 /// HDF5 file handle.
63 hid_t file = H5I_INVALID_HID;
64 /// The 'VTKHDF' root group within the file.
65 hid_t vtk = H5I_INVALID_HID;
66 /// Data transfer property list.
67 hid_t dxpl = H5P_DEFAULT;
68 /// The group to cell data (element attributes).
69 hid_t cell_data = H5I_INVALID_HID;
70 /// The group to store point data (e.g. grid functions).
71 hid_t point_data = H5I_INVALID_HID;
72
73 /// Compression level (-1 means disabled, 0 through 9 enabled). Default is 6.
74 int compression_level = 6;
75
76 /// Wrapper for storing dataset dimensions (max ndims is 2D in VTKHDF).
77 struct Dims
78 {
79 static constexpr size_t MAX_NDIMS = 2;
80 std::array<hsize_t, MAX_NDIMS> data = { }; // Zero initialized
81 int ndims = 0;
82 Dims() = default;
83 Dims(int ndims_) : ndims(ndims_) { MFEM_ASSERT(ndims <= MAX_NDIMS, ""); }
84 Dims(int ndims_, hsize_t val) : Dims(ndims_) { data.fill(val); }
85 template <typename T>
86 Dims(std::initializer_list<T> data_) : Dims(int(data_.size()))
87 { std::copy(data_.begin(), data_.end(), data.begin()); }
88 operator hsize_t*() { return data.data(); }
89 hsize_t &operator[](int i) { return data[i]; }
90 hsize_t TotalSize() const;
91 };
92
93 /// @name Needed for time-dependent data sets.
94 ///@{
95
96 /// The group to store time step information.
97 hid_t steps = H5I_INVALID_HID;
98
99 /// Number of time steps saved.
100 unsigned long nsteps = 0;
101
102 /// Keep track of the offsets into the data arrays at each time step.
103 struct Offsets
104 {
105 hsize_t current = 0;
106 hsize_t next = 0;
107 void Update(hsize_t offset)
108 {
109 current = next;
110 next += offset;
111 }
112 };
113
114 hsize_t part_offset; ///< Offset into the "NumberOf" arrays.
115 Offsets point_offsets; ///< Offsets into the point-sized arrays.
116 Offsets cell_offsets; ///< Offsets into the cell-sized arrays.
117 Offsets connectivity_offsets; ///< Offsets into the connectivity arrays.
118
119 /// Offsets into named point data arrays (for saved GridFunctions).
120 std::unordered_map<std::string,Offsets> point_data_offsets;
121
122 /// Track when the mesh has changed (enable reusing the previous saved mesh).
123 class MeshId
124 {
125 const Mesh *mesh_ptr = nullptr;
126 long sequence = -1;
127 long nodes_sequence = -1;
128 bool high_order = true;
129 int ref = -1;
130 public:
131 MeshId() = default;
132 /// @brief Is the given @a mesh different from the cached MeshId?
133 ///
134 /// The mesh is the same if the pointer is the same and sequence and nodes
135 /// sequence are both the same.
136 ///
137 /// If HasChanged() returns true, then new mesh data should be saved. If
138 /// it returns false, then the mesh doesn't need to be saved again.
139 bool HasChanged(const Mesh &mesh) const
140 {
141 if (mesh_ptr == &mesh)
142 {
143 return sequence != mesh.GetSequence() ||
144 nodes_sequence != mesh.GetNodesSequence();
145 }
146 return true;
147 }
148 /// @brief Is the given @a mesh different from the cached MeshId, taking
149 /// into account the order and refinement level of the previously saved
150 /// mesh.
151 bool HasChanged(const Mesh &mesh, bool high_order_, int ref_) const
152 {
153 return HasChanged(mesh) || high_order != high_order_ || ref != ref_;
154 }
155 /// Set the cached MeshId given @a mesh.
156 void Set(const Mesh &mesh, bool high_order_, int ref_)
157 {
158 mesh_ptr = &mesh;
159 sequence = mesh.GetSequence();
160 nodes_sequence = mesh.GetNodesSequence();
161 high_order = high_order_;
162 ref = ref_;
163 }
164 /// Return the refinement level of the previously saved mesh.
165 int GetRefinementLevel() const { return ref; }
166 };
167
168 /// The most recently saved MeshId.
169 MeshId mesh_id;
170
171 /// Number of points of most recently saved mesh.
172 hsize_t last_np = 0;
173
174 ///@}
175
176 /// Hold rank-offset and total size for parallel I/O.
177 struct OffsetTotal { size_t offset; size_t total; };
178
179 /// Ensure that the 'Steps' group and 'PointDataOffsets' subgroup exist.
180 void EnsureSteps();
181
182 /// @brief Ensure that the dataset named @a name in @a f exists.
183 ///
184 /// If the dataset does not exist, create it and return the ID. If the
185 /// dataset exists, open it and return the ID.
186 ///
187 /// The rank (number of dimensions) of the dataset is given by @a ndims and
188 /// its data type is given by @a type.
189 ///
190 /// If the dataset does not exist, it will initially have size @a dims.
191 /// Otherwise, it will be resized to append data of size @a dims, and @a dims
192 /// will be set to the new total size.
193 hid_t EnsureDataset(hid_t f, const std::string &name, hid_t type, Dims &dims);
194
195 /// @brief Ensure the named group is open, creating it if needed. Set @a
196 /// group to the ID.
197 ///
198 /// If @a group has already been opened, do nothing.
199 void EnsureGroup(const std::string &name, hid_t &group);
200
201 /// Appends data in parallel to the dataset named @a name in @a f.
202 ///
203 /// Data is appended along the zeroth dimension. Data of length @a locsize
204 /// will be written at offset @a offset. The sum over all MPI ranks of @a
205 /// locsize is equal to the zeroth dimension of @a globsize. The data is
206 /// written in row-major order.
207 template <typename T>
208 void AppendParData(hid_t f, const std::string &name, hsize_t locsize,
209 hsize_t offset, Dims globsize, T *data);
210
211 /// Appends data in parallel to the dataset named @a name in @a f.
212 ///
213 /// The input parameter @a dims is used to determine the rank of the data to
214 /// be written, and the extent of all but the first dimension. For example,
215 /// if dims = (0), then the written dimensions will be (N), where N is the
216 /// sum of the size of @a data across all MPI ranks. If dims = (0, m), then
217 /// the written dimensions will be (N/m, m). The data will be written in
218 /// row-major ordering.
219 ///
220 /// The row offset and total are returned.
221 template <typename T>
222 OffsetTotal AppendParVector(hid_t f, const std::string &name,
223 const std::vector<T> &data, Dims dims = Dims(1));
224
225 /// @brief Append a single value to the dataset named @a name in @a f.
226 template <typename T>
227 void AppendValue(const hid_t f, const std::string &name, T value);
228
229 /// Gather the value 'loc' from all MPI ranks and return the resulting array.
230 template <typename T>
231 std::vector<T> AllGather(const T loc) const;
232
233 /// @brief Returns the pair (offset, total) across MPI ranks given local data
234 /// 'loc'.
235 ///
236 /// If loc_i represents the value of loc on MPI rank i, then the returned
237 /// offset is the sum from i = 0 to R, where R is the current MPI rank. The
238 /// returned total is the sum over all MPI ranks.
239 ///
240 /// Requires performing 'gather all' operation.
241 OffsetTotal GetOffsetAndTotal(const size_t loc) const;
242
243 /// @brief Return true if the VTKHDF file is using MPI, false otherwise.
244 ///
245 /// The object is using MPI if the MPI communicator was passed to the
246 /// constructor. This is only possible if MFEM_USE_MPI is enabled. Even with
247 /// MPI enabled, a non-MPI VTKHDF object may be created.
248 bool UsingMpi() const;
249
250 /// Calls MPI_Barrier if in parallel, no-op otherwise.
251 void Barrier() const;
252
253 /// Return the HDF5 type ID corresponding to type @a T.
254 template <typename T> static hid_t GetTypeID();
255
256 /// Common setup (VTK group creation, etc.) for serial and parallel.
257 void SetupVTKHDF();
258
259 /// Create a new file (deleting existing file if needed).
260 void CreateFile(const std::string &filename, Restart restart);
261
262 /// Read the entire named dataset.
263 template <typename T>
264 std::vector<T> ReadDataset(const std::string &name) const;
265
266 /// Read a single value from the named dataset.
267 template <typename T>
268 T ReadValue(const std::string &name, hsize_t index) const;
269
270 /// Truncate the named dataset after position size in the first dimension.
271 void TruncateDataset(const std::string &name, hsize_t size);
272
273 /// Truncate all datasets on and after time @a t.
274 void Truncate(const real_t t);
275
276public:
277
278 /// @brief Create a new %VTKHDF file for serial I/O.
279 ///
280 /// If @a restart is enabled, then the file (if it exists) will be opened,
281 /// and time steps before the given time will be preserved. If @a restart
282 /// is not enabled, the file (if it exists) will be deleted, and a new file
283 /// will be created.
284 VTKHDF(const std::string &filename, Restart restart = Restart::Disabled());
285
286#ifdef MFEM_PARALLEL_HDF5
287 /// @brief Create a new %VTKHDF file for parallel I/O.
288 ///
289 /// If @a restart is enabled, then the file (if it exists) will be opened,
290 /// and time steps before the given time will be preserved. If @a restart
291 /// is not enabled, the file (if it exists) will be deleted, and a new file
292 /// will be created.
293 ///
294 /// If using restart mode, the file must have previously been saved with the
295 /// same number of MPI ranks.
296 VTKHDF(const std::string &filename, MPI_Comm comm_,
297 Restart restart = Restart::Disabled());
298#endif
299
300 /// @name Not copyable or movable.
301 ///@{
302 VTKHDF(const VTKHDF &) = delete;
303 VTKHDF(VTKHDF &&) = delete;
304 VTKHDF &operator=(const VTKHDF &) = delete;
305 VTKHDF &operator=(VTKHDF &&) = delete;
306 ///@}
307
308 /// Update the time step data after saving the mesh and grid functions.
309 void UpdateSteps(real_t t);
310
311 /// Disable zlib compression.
312 void DisableCompression() { compression_level = -1; }
313
314 /// @brief Enable zlib compression at the specified level.
315 ///
316 /// @a level must be between 0 and 9, in increasing order of compression.
317 void EnableCompression(int level = 6) { compression_level = level; }
318
319 /// @brief Save the mesh, appending as a new time step.
320 ///
321 /// If @a high_order is true, @a ref determines the polynomial degree of the
322 /// mesh elements (-1 indicates the same as the mesh nodes). If @a high_order
323 /// is false, the elements are uniformly subdivided according to @a ref.
324 template <typename FP_T = real_t>
325 void SaveMesh(const Mesh &mesh, bool high_order = true, int ref = -1);
326
327 /// Save the grid function with the given name, appending as a new time step.
328 template <typename FP_T = real_t>
329 void SaveGridFunction(const GridFunction &gf, const std::string &name);
330
331 /// Flush the file.
332 void Flush();
333
334 ///< Destructor. Close the file.
335 ~VTKHDF();
336};
337
338} // namespace mfem
339
340#endif // MFEM_USE_HDF5
341#endif // MFEM_VTKHDF
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:32
Mesh data type.
Definition mesh.hpp:65
Low-level class for writing VTKHDF data (for use in ParaView).
Definition vtkhdf.hpp:37
void DisableCompression()
Disable zlib compression.
Definition vtkhdf.hpp:312
VTKHDF(const VTKHDF &)=delete
void EnableCompression(int level=6)
Enable zlib compression at the specified level.
Definition vtkhdf.hpp:317
void SaveGridFunction(const GridFunction &gf, const std::string &name)
Save the grid function with the given name, appending as a new time step.
Definition vtkhdf.cpp:753
void SaveMesh(const Mesh &mesh, bool high_order=true, int ref=-1)
Save the mesh, appending as a new time step.
Definition vtkhdf.cpp:534
VTKHDF & operator=(VTKHDF &&)=delete
void UpdateSteps(real_t t)
Update the time step data after saving the mesh and grid functions.
Definition vtkhdf.cpp:507
VTKHDF(const std::string &filename, Restart restart=Restart::Disabled())
Create a new VTKHDF file for serial I/O.
Definition vtkhdf.cpp:460
VTKHDF & operator=(const VTKHDF &)=delete
VTKHDF(VTKHDF &&)=delete
void Flush()
Flush the file.
Definition vtkhdf.cpp:789
int index(int i, int j, int nx, int ny)
Definition life.cpp:236
mfem::real_t real_t
float real_t
Definition config.hpp:46
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
Helper used in VTKHDF::VTKHDF for enabling/disabling restart mode.
Definition vtkhdf.hpp:41
Restart(real_t time_)
Definition vtkhdf.hpp:45
Restart(bool enabled_, real_t time_)
Definition vtkhdf.hpp:46
static Restart Disabled()
Definition vtkhdf.hpp:47