DNDSR 0.2.1
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
SerializerH5.cpp
Go to the documentation of this file.
1#include "SerializerH5.hpp"
2#include <algorithm>
3#include <string>
4#include <fmt/core.h>
6
7namespace DNDS::Serializer
8{
9#define H5CHECK_Set DNDS_assert_info(herr >= 0, "H5 setting err")
10#define H5CHECK_Close DNDS_assert_info(herr >= 0, "H5 closing err")
11#define H5CHECK_Iter DNDS_assert_info(herr >= 0, "H5 iteration err")
12 /**
13 * @brief returns good path and if the path is absolute
14 */
15 static bool processPath(std::vector<std::string> &pth)
16 {
17 bool ifAbs = pth.empty() || pth[0].empty();
18 pth.erase(
19 std::remove_if(pth.begin(), pth.end(), [](const std::string &v)
20 { return v.empty(); }), // only keep non-zero sized path names
21 pth.end());
22
23 return ifAbs;
24 }
25
26 static std::string constructPath(std::vector<std::string> &pth)
27 {
28 std::string ret;
29 for (auto &name : pth)
30 ret.append(std::string("/") + name);
31 return ret;
32 }
33
34 /**
35 * @brief Get the Group Of File If Exist
36 *
37 * @param file_id
38 * @param read
39 * @param groupName
40 * @return hid_t group_id, need releasing
41 * @warning remember releasing returned group_id!
42 */
43 static hid_t GetGroupOfFileIfExist(hid_t file_id, bool read, const std::string &groupName, bool coll_on_meta)
44 {
45 hid_t group_id{-1};
46 herr_t herr{0};
47
48 auto pth = splitSString(groupName, '/');
49 bool isAbs = processPath(pth);
50 // std::cout << groupName << std::endl;
51 DNDS_assert_info(isAbs, fmt::format("groupName: {}", groupName));
52 if (pth.empty()) // is root itself
53 {
54 group_id = H5Gopen(file_id, "/", H5P_DEFAULT);
55 DNDS_assert(group_id >= 0);
56 }
57 for (int i = 0; i < pth.size(); i++)
58 {
59 std::vector<std::string> pth_parent;
60 for (int j = 0; j <= i; j++)
61 pth_parent.push_back(pth[j]);
62 auto parentGroupName = constructPath(pth_parent);
63 hid_t lapl_id = H5Pcreate(H5P_LINK_ACCESS);
64 DNDS_assert(lapl_id >= 0);
65 if (coll_on_meta)
66 herr = H5Pset_all_coll_metadata_ops(lapl_id, true), H5CHECK_Set;
67 htri_t group_exists = H5Lexists(file_id, parentGroupName.c_str(), lapl_id);
68 herr = H5Pclose(lapl_id), H5CHECK_Close;
69 if (group_exists > 0)
70 {
71 // If the group exists, open it if it the deepest level
72 if (i == pth.size() - 1)
73 {
74 hid_t gapl_id = H5Pcreate(H5P_GROUP_ACCESS);
75 DNDS_assert(gapl_id >= 0);
76 if (coll_on_meta)
77 herr = H5Pset_all_coll_metadata_ops(gapl_id, true), H5CHECK_Set;
78 group_id = H5Gopen(file_id, parentGroupName.c_str(), gapl_id);
79 herr = H5Pclose(gapl_id), H5CHECK_Close;
80 // std::cout << parentGroupName << " exists" << std::endl;
81 }
82 else
83 group_id = -1;
84 }
85 else
86 {
87 DNDS_assert_info(!read, "file is read only, cannot create new group: " + groupName);
88
89 // If the group does not exist, create it
90 hid_t gapl_id = H5Pcreate(H5P_GROUP_ACCESS);
91 DNDS_assert(gapl_id >= 0);
92 if (coll_on_meta)
93 herr = H5Pset_all_coll_metadata_ops(gapl_id, true), H5CHECK_Set;
94 group_id = H5Gcreate(file_id, parentGroupName.c_str(), H5P_DEFAULT, H5P_DEFAULT, gapl_id);
95 DNDS_assert(group_id >= 0);
96 herr = H5Pclose(gapl_id), H5CHECK_Close;
97 // std::cout << parentGroupName << " created" << std::endl;
98 if (i < pth.size() - 1)
99 {
100 herr = H5Gclose(group_id), H5CHECK_Close;
101 group_id = -1;
102 }
103 }
104 }
105
106 return group_id;
107 }
108
109 void SerializerH5::OpenFile(const std::string &fName, bool read)
110 {
111 reading = read;
112 cP = "";
113 herr_t herr{0};
114 hid_t plist_id = H5Pcreate(H5P_FILE_ACCESS);
115 herr = H5Pset_fapl_mpio(plist_id, commDup, MPI_INFO_NULL), H5CHECK_Set; // Set up file access property list with parallel I/O access
116 herr = H5Pset_all_coll_metadata_ops(plist_id, true), H5CHECK_Set;
117 herr = H5Pset_coll_metadata_write(plist_id, true), H5CHECK_Set;
118 if (read)
119 h5file = H5Fopen(fName.c_str(), H5F_ACC_RDONLY, plist_id);
120 else
121 h5file = H5Fcreate(fName.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, plist_id);
122 DNDS_assert_info(H5I_INVALID_HID != h5file, fmt::format(" attempted to {} file [{}]", read ? "read" : "write", fName));
123 herr = H5Pclose(plist_id), H5CHECK_Close;
124 }
130 {
131 cPathSplit.clear();
132 dedupClear();
133 cP.clear();
134
135 if (H5I_INVALID_HID != h5file)
136 H5Fclose(h5file), h5file = H5I_INVALID_HID;
137 }
138 void SerializerH5::CreatePath(const std::string &p)
139 {
140 herr_t herr{0};
141 auto pth = splitSString(p, '/');
142 bool isAbs = processPath(pth);
143 std::vector<std::string> newPath = cPathSplit;
144 if (isAbs)
145 newPath = std::move(pth);
146 else
147 for (auto &name : pth)
148 newPath.push_back(name);
149 std::string nP = constructPath(cPathSplit);
150 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, nP, collectiveMetadataRW);
151 herr = H5Gclose(group_id), H5CHECK_Close;
152 }
153 void SerializerH5::GoToPath(const std::string &p)
154 {
155 auto pth = splitSString(p, '/');
156 bool isAbs = processPath(pth);
157 if (isAbs)
158 cPathSplit = std::move(pth);
159 else
160 for (auto &name : pth)
161 cPathSplit.push_back(name);
162 cP = constructPath(cPathSplit);
163 // DNDS_assert(jObj[nlohmann::json::json_pointer(cP)].is_object());
164 }
166 {
167 return cP;
168 }
169
170 // Structure to hold collected HDF5 information
172 {
173 std::vector<std::string> groups;
174 std::vector<std::string> datasets;
175 std::vector<std::string> attributes;
176 };
177
179 {
180 // TraverseData is a short-lived per-call aggregate passed through the
181 // HDF5 H5Literate callback; the reference is intentional so that the
182 // callback mutates the caller's H5Contents directly.
183 // NOLINTNEXTLINE(cppcoreguidelines-avoid-const-or-ref-data-members)
185 std::string current_path;
187 [[nodiscard]] std::string get_indent() const
188 {
189 // Brace-init `{n, ' '}` is ambiguous with std::string's
190 // initializer_list<char> ctor and triggers -Wnarrowing; keep
191 // the explicit std::string(n, ch) form.
192 // NOLINTNEXTLINE(modernize-return-braced-init-list)
193 return std::string(std::count(current_path.begin(), current_path.end(), '/') * 2, ' ');
194 }
195 };
196
197 static herr_t link_iterate_cb(hid_t group_id, const char *name, const H5L_info_t *info, void *op_data)
198 {
199 herr_t herr = 0;
200 auto *data = static_cast<TraverseData *>(op_data);
201 bool coll_on_meta = data->coll_on_meta;
202 std::string full_name = data->current_path + "/" + name;
203 // We only care about hard links that represent HDF5 objects
204 // std::cout << "name is " << name << std::endl;
205 if (info->type == H5L_TYPE_HARD)
206 {
207 data->contents.groups.emplace_back(name);
208 return 0;
209 //! now it seems obj_info is not correctly retrieved with type=Unknown
210 //! TODO: fix this and get actual object type
211 H5O_info_t obj_info;
212 hid_t lapl_id = H5Pcreate(H5P_LINK_ACCESS);
213 DNDS_assert(lapl_id >= 0);
214 herr = H5Pset_all_coll_metadata_ops(lapl_id, true), H5CHECK_Set;
215 // if (coll_on_meta)
216 // herr = H5Pset_all_coll_metadata_ops(lapl_id, true), H5CHECK_Set;
217 if (H5Oget_info_by_name(group_id, name, &obj_info, H5P_DEFAULT, lapl_id) >= 0)
218 {
219 log() << "name is " << name << ", \n"
220 << obj_info.type << "\n"
221 << obj_info.mtime << "\n"
222 << obj_info.rc << "\n"
223 << obj_info.fileno << "\n"
224 << std::endl;
225 switch (obj_info.type)
226 {
227 case H5O_TYPE_GROUP:
228 {
229 data->contents.groups.emplace_back(name);
230 // std::cout << data->get_indent() << "Group: " << full_name << std::endl;
231 break;
232 }
233 case H5O_TYPE_DATASET:
234 {
235 data->contents.datasets.emplace_back(name);
236 // std::cout << data->get_indent() << "Dataset: " << full_name << std::endl;
237 break;
238 }
239 default:
240 // Other types (e.g., symbols, but less common for direct user interaction)
241 break;
242 }
243 }
244 else
245 {
246 log() << "Error getting object info for: " << full_name << std::endl;
247 H5Eprint2(H5E_DEFAULT, stderr);
248 }
249 herr = H5Pclose(lapl_id), H5CHECK_Close;
250 }
251 else
252 {
253 // Handle soft links or external links if their names are also desired
254 // std::cout << data->get_indent() << "Link: " << full_name << " (Non-hard link)" << std::endl;
255 }
256 return 0; // Continue iteration
257 }
258
259 // Callback for H5Aiterate (iterating attributes)
260 static herr_t attribute_iterate_cb(hid_t obj_id, const char *attr_name, const H5A_info_t *info, void *op_data)
261 {
262 auto *data = static_cast<TraverseData *>(op_data);
263 data->contents.attributes.emplace_back(attr_name);
264 // std::string full_attr_name = data->current_path + "@" + attr_name; // Common convention for attribute paths
265 // std::cout << data->get_indent() << " Attribute: " << full_attr_name << std::endl; // Indent more for attributes
266 return 0; // Continue iteration
267 }
268
269 std::set<std::string> SerializerH5::ListCurrentPath()
270 {
271 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, collectiveMetadataRW);
272 H5Contents contents;
273 TraverseData data{contents, cP, collectiveMetadataRW};
274 herr_t herr{0};
275 herr = H5Aiterate(group_id, H5_INDEX_NAME, H5_ITER_INC, nullptr, attribute_iterate_cb, &data), H5CHECK_Iter;
276 herr = H5Literate(group_id, H5_INDEX_NAME, H5_ITER_INC, nullptr, link_iterate_cb, &data), H5CHECK_Iter;
277 H5Gclose(group_id), H5CHECK_Close; // don't forget this
278 std::set<std::string> ret;
279 for (auto &v : contents.attributes)
280 ret.insert(v);
281 for (auto &v : contents.groups)
282 ret.insert(v);
283 for (auto &v : contents.datasets)
284 ret.insert(v);
285 return ret;
286 }
287
288 template <typename T>
289 // using T = int;
290 void WriteAttributeScalar(const std::string &name, const T &v,
291 hid_t h5file, bool reading, const std::string &cP,
292 bool coll_on_meta)
293 {
294 hid_t T_H5TYPE = -1;
295 if constexpr (std::is_same_v<T, int>)
296 T_H5TYPE = H5T_NATIVE_INT;
297 else if constexpr (std::is_same_v<T, index>)
298 T_H5TYPE = DNDS_H5T_INDEX();
299 else if constexpr (std::is_same_v<T, rowsize>)
300 T_H5TYPE = DNDS_H5T_ROWSIZE();
301 else if constexpr (std::is_same_v<T, real>)
302 T_H5TYPE = DNDS_H5T_REAL();
303 else if constexpr (std::is_same_v<T, std::string>)
304 ; // pass
305 else
306 static_assert(std::is_same_v<T, real>);
307
308 // `vV` is addressed via `&vV` in the non-string `if constexpr` branch
309 // below; clang-tidy's check misses the template instantiation for
310 // `T != std::string`.
311 // NOLINTNEXTLINE(performance-unnecessary-copy-initialization)
312 T vV = v;
313
314 if constexpr (!std::is_same_v<T, std::string>)
315 {
316 herr_t herr{0};
317 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, coll_on_meta);
318 hsize_t attr_size = 1;
319 hid_t attr_space = H5Screate(H5S_SCALAR);
320 hid_t aapl_id = H5Pcreate(H5P_ATTRIBUTE_ACCESS);
321 if (coll_on_meta)
322 herr = H5Pset_all_coll_metadata_ops(aapl_id, true), H5CHECK_Set;
323 hid_t attr_id = H5Acreate(group_id, name.c_str(), T_H5TYPE, attr_space, H5P_DEFAULT, aapl_id);
324 herr = H5Awrite(attr_id, T_H5TYPE, &vV), H5CHECK_Set;
325 H5Aclose(attr_id);
326 H5Pclose(aapl_id);
327 H5Sclose(attr_space);
328 H5Gclose(group_id);
329 }
330 else
331 {
332 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, coll_on_meta);
333
334 herr_t herr{0};
335 hid_t aapl_id = H5Pcreate(H5P_ATTRIBUTE_ACCESS);
336 if (coll_on_meta)
337 herr = H5Pset_all_coll_metadata_ops(aapl_id, true), H5CHECK_Set;
338 hid_t scalar_space = H5Screate(H5S_SCALAR);
339 hid_t string_type = H5Tcreate(H5T_STRING, v.length());
340 herr = H5Tset_strpad(string_type, H5T_STR_NULLPAD), H5CHECK_Set;
341 hid_t type_attr_id = H5Acreate(group_id, name.c_str(), string_type, scalar_space, H5P_DEFAULT, aapl_id);
342 DNDS_assert_info(type_attr_id >= 0, fmt::format("{}, {}, {}", name, v, cP));
343 herr = H5Awrite(type_attr_id, string_type, v.data()), H5CHECK_Set;
344
345 herr = H5Aclose(type_attr_id), H5CHECK_Close;
346 herr = H5Tclose(string_type), H5CHECK_Close;
347 herr = H5Pclose(aapl_id), H5CHECK_Close;
348 herr = H5Sclose(scalar_space), H5CHECK_Close;
349 herr = H5Gclose(group_id), H5CHECK_Close;
350 }
351 }
352
353 void SerializerH5::WriteInt(const std::string &name, int v)
354 {
355 WriteAttributeScalar<int>(name, v, h5file, reading, cP, collectiveMetadataRW);
356 }
357 void SerializerH5::WriteIndex(const std::string &name, index v)
358 {
359 WriteAttributeScalar<index>(name, v, h5file, reading, cP, collectiveMetadataRW);
360 }
361 void SerializerH5::WriteReal(const std::string &name, real v)
362 {
363 WriteAttributeScalar<real>(name, v, h5file, reading, cP, collectiveMetadataRW);
364 }
365
366 static void H5_WriteDataset(hid_t loc, const char *name, index nGlobal, index nOffset, index nLocal,
367 hid_t file_dataType, hid_t mem_dataType, hid_t dxpl_id, hid_t dapl_id, int64_t chunksize, int deflateLevel,
368 const void *buf, int dim2 = -1)
369 {
370 int herr{0};
371 DNDS_assert_info(nGlobal >= 0 && nLocal >= 0 && nOffset >= 0,
372 fmt::format("{},{},{}", nGlobal, nLocal, nOffset));
373 int rank = dim2 >= 0 ? 2 : 1;
374 std::array<hsize_t, 2> ranksFull{hsize_t(nGlobal), hsize_t(dim2)};
375 std::array<hsize_t, 2> ranksFullUnlim{chunksize > 0 ? H5S_UNLIMITED : hsize_t(nGlobal), hsize_t(dim2)};
376 std::array<hsize_t, 2> offset{hsize_t(nOffset), 0};
377 std::array<hsize_t, 2> siz{hsize_t(nLocal), hsize_t(dim2)};
378 hid_t memSpace = H5Screate_simple(rank, siz.data(), nullptr);
379 hid_t fileSpace = H5Screate_simple(rank, ranksFull.data(), ranksFullUnlim.data());
380 std::array<hsize_t, 2> chunk_dims{hsize_t(chunksize > 0 ? chunksize : 0), dim2 >= 0 ? hsize_t(dim2) : 0};
381 hid_t dcpl_id = H5Pcreate(H5P_DATASET_CREATE);
382 if (chunk_dims[0] > 0)
383 herr = H5Pset_chunk(dcpl_id, dim2 > 0 ? 2 : 1, chunk_dims.data()), H5CHECK_Set;
384#ifdef H5_HAVE_FILTER_DEFLATE
385 if (deflateLevel > 0)
386 herr = H5Pset_deflate(dcpl_id, deflateLevel), H5CHECK_Set;
387#endif
388
389 hid_t dset_id = H5Dcreate(loc, name, file_dataType, fileSpace, H5P_DEFAULT, dcpl_id, dapl_id);
390 DNDS_assert_info(H5I_INVALID_HID != dset_id, "dataset create failed");
391 herr = H5Sclose(fileSpace);
392 fileSpace = H5Dget_space(dset_id);
393 herr |= H5Sselect_hyperslab(fileSpace, H5S_SELECT_SET, offset.data(), nullptr, siz.data(), nullptr);
394 herr |= H5Dwrite(dset_id, mem_dataType, memSpace, fileSpace, dxpl_id, buf);
395 herr |= H5Dclose(dset_id);
396 herr |= H5Pclose(dcpl_id);
397 herr |= H5Sclose(fileSpace);
398 herr |= H5Sclose(memSpace);
399 DNDS_assert_info(herr >= 0,
400 "h5 error " + fmt::format(
401 "nGlobal {}, nOffset {}, nLocal {}, name {}",
402 nGlobal, nOffset, nLocal, name));
403 } //! currently a direct copy
404
405 template <typename T = index>
406 // using T = index;
407 static void WriteDataVector(const std::string &name, const T *v, size_t size, ArrayGlobalOffset offset, int64_t chunksize, int deflateLevel,
408 hid_t h5file, bool reading, const std::string &cP, const MPIInfo &mpi, MPI_Comm commDup,
409 bool coll_on_meta, bool coll_on_data)
410 {
411 hid_t T_H5TYPE = H5T_NATIVE_INT;
412 if constexpr (std::is_same_v<T, index>)
413 T_H5TYPE = DNDS_H5T_INDEX();
414 else if constexpr (std::is_same_v<T, rowsize>)
415 T_H5TYPE = DNDS_H5T_ROWSIZE();
416 else if constexpr (std::is_same_v<T, real>)
417 T_H5TYPE = DNDS_H5T_REAL();
418 else if constexpr (std::is_same_v<T, uint8_t>)
419 T_H5TYPE = H5T_NATIVE_UINT8;
420 else
421 static_assert(std::is_same_v<T, uint8_t>);
422
423 herr_t herr{0};
424 hid_t dxpl_id = H5Pcreate(H5P_DATASET_XFER);
425 //! is this necessary? seems even global-unique vector needs this if deflate is on!
426 if (((offset.isDist() || offset == ArrayGlobalOffset_Parts) && (coll_on_data)) || deflateLevel > 0)
427 herr = H5Pset_dxpl_mpio(dxpl_id, H5FD_MPIO_COLLECTIVE), H5CHECK_Set;
428 else
429 herr = H5Pset_dxpl_mpio(dxpl_id, H5FD_MPIO_INDEPENDENT), H5CHECK_Set;
430 DNDS_assert_info(herr >= 0, "h5 error");
431 // if is dist array, we use coll metadata
432 hid_t dapl_id = H5Pcreate(H5P_DATASET_ACCESS);
433 if (coll_on_meta)
434 herr = H5Pset_all_coll_metadata_ops(dapl_id, true), H5CHECK_Set;
435
436 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, coll_on_meta);
437 if (offset == ArrayGlobalOffset_One)
438 H5_WriteDataset(group_id, name.c_str(), size, 0, mpi.rank == 0 ? size : 0,
439 T_H5TYPE, T_H5TYPE, dxpl_id, dapl_id, chunksize, deflateLevel,
440 v);
441 else if (offset == ArrayGlobalOffset_Parts) // now we force it to be
442 {
443 uint64_t sizeU{size}, sizeOff{0}, sizeGlobal{0};
444 MPI::Scan(&sizeU, &sizeOff, 1, MPI_UINT64_T, MPI_SUM, commDup);
445 sizeGlobal = sizeOff;
446 MPI::Bcast(&sizeGlobal, 1, MPI_UINT64_T, mpi.size - 1, commDup);
447 sizeOff -= sizeU;
448
449 H5_WriteDataset(group_id, name.c_str(), sizeGlobal, sizeOff, size,
450 T_H5TYPE, T_H5TYPE, dxpl_id, dapl_id, chunksize, deflateLevel,
451 v);
452 // rank offset array
453 std::array<index, 2> offsetC = {index(sizeOff), index(sizeGlobal)};
454 H5_WriteDataset(group_id, (name + "::rank_offsets").c_str(), mpi.size + 1, mpi.rank, (mpi.rank == mpi.size - 1) ? 2 : 1,
455 DNDS_H5T_INDEX(), DNDS_H5T_INDEX(), dxpl_id, dapl_id, chunksize, deflateLevel,
456 offsetC.data());
457 }
458 else if (offset.isDist())
459 {
460
461 H5_WriteDataset(group_id, name.c_str(), offset.size(), offset.offset(), size,
462 T_H5TYPE, T_H5TYPE, dxpl_id, dapl_id, chunksize, deflateLevel,
463 v);
464 // rank offset array
465 std::array<index, 2> offsetC = {offset.offset(), offset.size()};
466 H5_WriteDataset(group_id, (name + "::rank_offsets").c_str(), mpi.size + 1, mpi.rank, (mpi.rank == mpi.size - 1) ? 2 : 1,
467 DNDS_H5T_INDEX(), DNDS_H5T_INDEX(), dxpl_id, dapl_id, chunksize, deflateLevel,
468 offsetC.data());
469 }
470 else
471 DNDS_assert_info(false, "offset ill-formed");
472
473 herr = H5Pclose(dxpl_id), H5CHECK_Close;
474 herr = H5Pclose(dapl_id), H5CHECK_Close;
475 herr = H5Gclose(group_id), H5CHECK_Close;
476 }
477
478 void SerializerH5::WriteIndexVector(const std::string &name, const std::vector<index> &v, ArrayGlobalOffset offset)
479 {
480 WriteDataVector<index>(name, v.data(), v.size(), offset, chunksize, deflateLevel, h5file, reading, cP, mpi, commDup, collectiveMetadataRW, collectiveDataRW);
481 }
482 void SerializerH5::WriteRowsizeVector(const std::string &name, const std::vector<rowsize> &v, ArrayGlobalOffset offset)
483 {
484 WriteDataVector<rowsize>(name, v.data(), v.size(), offset, chunksize, deflateLevel, h5file, reading, cP, mpi, commDup, collectiveMetadataRW, collectiveDataRW);
485 }
486 void SerializerH5::WriteRealVector(const std::string &name, const std::vector<real> &v, ArrayGlobalOffset offset)
487 {
488 WriteDataVector<real>(name, v.data(), v.size(), offset, chunksize, deflateLevel, h5file, reading, cP, mpi, commDup, collectiveMetadataRW, collectiveDataRW);
489 }
490 void SerializerH5::WriteString(const std::string &name, const std::string &v)
491 {
492 WriteAttributeScalar<std::string>(name, v, h5file, reading, cP, collectiveMetadataRW);
493 }
495 {
496 std::string refPath;
497 index CanNotShare = dedupLookup(v, refPath) ? 0 : 1;
498 MPI::AllreduceOneIndex(CanNotShare, MPI_MAX, MPIInfo{commDup, mpi.rank, mpi.size});
499 if (!CanNotShare)
500 this->WriteString(name + "::ref", refPath);
501 else
502 {
503 WriteDataVector<index>(name, v->data(), v->size(), offset, chunksize, deflateLevel, h5file, reading, cP, mpi, commDup, collectiveMetadataRW, collectiveDataRW);
504 dedupRegister(v, cP + "/" + name);
505 }
506 }
508 {
509 std::string refPath;
510 index CanNotShare = dedupLookup(v, refPath) ? 0 : 1;
511 MPI::AllreduceOneIndex(CanNotShare, MPI_MAX, MPIInfo{commDup, mpi.rank, mpi.size});
512 if (!CanNotShare)
513 this->WriteString(name + "::ref", refPath);
514 else
515 {
516 WriteDataVector<rowsize>(name, v->data(), v->size(), offset, chunksize, deflateLevel, h5file, reading, cP, mpi, commDup, collectiveMetadataRW, collectiveDataRW);
517 dedupRegister(v, cP + "/" + name);
518 }
519 }
520
521 template <typename T>
522 // using T = int;
523 void ReadAttributeScalar(const std::string &name, T &v,
524 hid_t h5file, bool reading, const std::string &cP,
525 bool coll_on_meta)
526 {
527 hid_t T_H5TYPE = -1;
528 if constexpr (std::is_same_v<T, int>)
529 T_H5TYPE = H5T_NATIVE_INT;
530 else if constexpr (std::is_same_v<T, index>)
531 T_H5TYPE = DNDS_H5T_INDEX();
532 else if constexpr (std::is_same_v<T, rowsize>)
533 T_H5TYPE = DNDS_H5T_ROWSIZE();
534 else if constexpr (std::is_same_v<T, real>)
535 T_H5TYPE = DNDS_H5T_REAL();
536 else if constexpr (std::is_same_v<T, std::string>)
537 ; // pass
538 else
539 static_assert(std::is_same_v<T, real>);
540
541 // `vV` is addressed via `&vV` in the non-string `if constexpr` branch
542 // below; clang-tidy's check misses the template instantiation for
543 // `T != std::string`.
544 // NOLINTNEXTLINE(performance-unnecessary-copy-initialization)
545 T vV = v;
546
547 if constexpr (!std::is_same_v<T, std::string>)
548 {
549 herr_t herr{0};
550 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, coll_on_meta);
551 hid_t aapl_id = H5Pcreate(H5P_ATTRIBUTE_ACCESS);
552 herr = H5Pset_all_coll_metadata_ops(aapl_id, true), H5CHECK_Set;
553 // herr = H5Pset_coll_metadata_write(aapl_id, true), H5CHECK_Set; // setting to present in dapl
554 hid_t attr_id = H5Aopen(group_id, name.c_str(), aapl_id);
555 DNDS_assert_info(attr_id >= 0, fmt::format("attempting to open attribute [{}/{}] failed", cP, name));
556 hid_t attr_space = H5Aget_space(attr_id);
557 int ndims = H5Sget_simple_extent_ndims(attr_space);
558 DNDS_assert_info(ndims == 0, fmt::format("attempting to read attribute [{}/{}] which is not scalar", cP, name));
559 herr = H5Aread(attr_id, T_H5TYPE, &v), H5CHECK_Set;
560 herr = H5Aclose(attr_id), H5CHECK_Close;
561 herr = H5Pclose(aapl_id), H5CHECK_Close;
562 herr = H5Sclose(attr_space), H5CHECK_Close;
563 herr = H5Gclose(group_id), H5CHECK_Close;
564 }
565 else
566 {
567 herr_t herr{0};
568 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, coll_on_meta);
569 hid_t aapl_id = H5Pcreate(H5P_ATTRIBUTE_ACCESS);
570 if (coll_on_meta)
571 {
572 herr = H5Pset_all_coll_metadata_ops(aapl_id, true), H5CHECK_Set;
573 // herr = H5Pset_coll_metadata_write(aapl_id, true), H5CHECK_Set; // setting not present in aapl
574 }
575 hid_t attr_id = H5Aopen(group_id, name.c_str(), aapl_id);
576 DNDS_assert_info(attr_id >= 0, fmt::format("attempting to open attribute [{}/{}] failed", cP, name));
577 hid_t attr_space = H5Aget_space(attr_id);
578 int ndims = H5Sget_simple_extent_ndims(attr_space);
579 DNDS_assert_info(ndims == 0, fmt::format("attempting to read attribute [{}/{}] which is not scalar", cP, name));
580 hid_t dtype_id = H5Aget_type(attr_id);
581 bool is_varlen = H5Tis_variable_str(dtype_id);
582
583 if (is_varlen)
584 {
585 // Variable-length string: HDF5 will allocate memory
586 char *attr_value = nullptr;
587 // H5Aread wants `void *buf`; the multi-level-implicit-pointer
588 // check requires an explicit reinterpret to silence the
589 // `char ** -> void *` conversion inside the argument parens.
590 // NOLINTNEXTLINE(bugprone-multi-level-implicit-pointer-conversion)
591 herr = H5Aread(attr_id, dtype_id, &attr_value), H5CHECK_Set;
592 // std::cout << "Read Attribute (Variable-Length): " << attr_value << "\n";
593 v = attr_value; // copy as null-terminated string
594 H5free_memory(attr_value); // Must free memory manually
595 }
596 else
597 {
598 // Fixed-length string: Allocate buffer based on datatype size
599 size_t size = H5Tget_size(dtype_id);
600 std::vector<char> buffer(size + 1, '\0'); // +1 for null terminator
601 herr = H5Aread(attr_id, dtype_id, buffer.data()), H5CHECK_Set;
602 // std::cout << "Read Attribute (Fixed-Length): " << buffer.data() << "\n";
603 v = buffer.data(); // copy as null-terminated string
604 }
605
606 herr = H5Tclose(dtype_id), H5CHECK_Close;
607 herr = H5Aclose(attr_id), H5CHECK_Close;
608 herr = H5Pclose(aapl_id), H5CHECK_Close;
609 herr = H5Sclose(attr_space), H5CHECK_Close;
610 herr = H5Gclose(group_id), H5CHECK_Close;
611 }
612 }
613
614 void SerializerH5::ReadInt(const std::string &name, int &v)
615 {
616 ReadAttributeScalar<int>(name, v, h5file, reading, cP, collectiveMetadataRW);
617 }
618 void SerializerH5::ReadIndex(const std::string &name, index &v)
619 {
620 ReadAttributeScalar<index>(name, v, h5file, reading, cP, collectiveMetadataRW);
621 }
622 void SerializerH5::ReadReal(const std::string &name, real &v)
623 {
624 ReadAttributeScalar<real>(name, v, h5file, reading, cP, collectiveMetadataRW);
625 }
626
627 /**
628 * @brief
629 *
630 * when #buf == 2, only set nGlobal and dim2
631 *
632 * @param loc
633 * @param name
634 * @param nGlobal
635 * @param nOffset
636 * @param nLocal
637 * @param file_dataType
638 * @param mem_dataType
639 * @param plist_id
640 * @param dcpl_id
641 * @param dapl_id
642 * @param buf
643 * @param dim2
644 */
645 // When dxpl_id uses collective I/O (H5FD_MPIO_COLLECTIVE), every rank MUST
646 // call H5Dread even if nLocal==0. Callers must provide a valid (possibly
647 // dummy) buf pointer; do NOT skip this function based on nLocal==0.
648 static void H5_ReadDataset(hid_t loc, const char *name, index &nGlobal, index nOffset, index nLocal,
649 hid_t mem_dataType, hid_t dxpl_id, hid_t dapl_id,
650 void *buf, int &dim2)
651 {
652 int herr{0};
653 hid_t dset_id = H5Dopen(loc, name, dapl_id);
654 DNDS_assert_info(dset_id >= 0, fmt::format("dataset [{}] open failed", name));
655 hid_t fileSpace = H5Dget_space(dset_id);
656 DNDS_assert_info(fileSpace >= 0, fmt::format("dataset [{}] filespace open failed", name));
657 int ndims = H5Sget_simple_extent_ndims(fileSpace);
658 DNDS_assert_info(ndims == 1 || ndims == 2, fmt::format("dataset [{}] not having 1 or 2 dims!", name));
659 std::array<hsize_t, 2> sizes{};
660 ndims = H5Sget_simple_extent_dims(fileSpace, sizes.data(), nullptr);
661 if (ndims == 2)
662 dim2 = sizes[1];
663 else
664 dim2 = -1;
665 nGlobal = sizes[0];
666
667 if (buf != nullptr)
668 {
669 int rank = dim2 >= 0 ? 2 : 1;
670 std::array<hsize_t, 2> offset{hsize_t(nOffset), 0};
671 std::array<hsize_t, 2> siz{hsize_t(nLocal), hsize_t(dim2)};
672 hid_t memSpace = H5Screate_simple(rank, siz.data(), nullptr);
673 DNDS_assert(memSpace > 0);
674 herr = H5Sselect_hyperslab(fileSpace, H5S_SELECT_SET, offset.data(), nullptr, siz.data(), nullptr), H5CHECK_Set;
675 herr = H5Dread(dset_id, mem_dataType, memSpace, fileSpace, dxpl_id, buf), H5CHECK_Set;
676 herr = H5Sclose(memSpace), H5CHECK_Close;
677 }
678 herr = H5Dclose(dset_id), H5CHECK_Close;
679 herr = H5Sclose(fileSpace), H5CHECK_Close;
680 }
681
682 template <typename T = index>
683 // using T = index;
684 static void ReadDataVector(const std::string &name, T *v, size_t &size, ArrayGlobalOffset &offset,
685 hid_t h5file, bool reading, const std::string &cP, const MPIInfo &mpi,
686 bool coll_on_meta, bool coll_on_data)
687 {
688 hid_t T_H5TYPE = H5T_NATIVE_INT;
689 if constexpr (std::is_same_v<T, index>)
690 T_H5TYPE = DNDS_H5T_INDEX();
691 else if constexpr (std::is_same_v<T, rowsize>)
692 T_H5TYPE = DNDS_H5T_ROWSIZE();
693 else if constexpr (std::is_same_v<T, real>)
694 T_H5TYPE = DNDS_H5T_REAL();
695 else if constexpr (std::is_same_v<T, uint8_t>)
696 T_H5TYPE = H5T_NATIVE_UINT8;
697 else
698 static_assert(std::is_same_v<T, uint8_t>);
699
700 herr_t herr{0};
701 hid_t dxpl_id = H5Pcreate(H5P_DATASET_XFER);
702 if ((offset.isDist() || offset == ArrayGlobalOffset_One) && coll_on_data) //! is this necessary?
703 herr = H5Pset_dxpl_mpio(dxpl_id, H5FD_MPIO_COLLECTIVE), H5CHECK_Set;
704 else
705 herr = H5Pset_dxpl_mpio(dxpl_id, H5FD_MPIO_INDEPENDENT), H5CHECK_Set;
706
707 // if is dist array, we use coll metadata
708 hid_t dapl_id = H5Pcreate(H5P_DATASET_ACCESS);
709 if (coll_on_meta)
710 {
711 herr = H5Pset_all_coll_metadata_ops(dapl_id, offset.isDist() || offset == ArrayGlobalOffset_One), H5CHECK_Set;
712 // herr = H5Pset_coll_metadata_write(dapl_id, offset.isDist() || offset == ArrayGlobalOffset_One), H5CHECK_Set; // tsetting to present in dapl
713 }
714 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, coll_on_meta);
715
716 if (offset == ArrayGlobalOffset_Unknown) // need detection, auto ArrayGlobalOffset_Parts or distributed array
717 {
718 DNDS_assert(v == nullptr); // must be a size-query call
719
720 hid_t lapl_id = H5Pcreate(H5P_LINK_ACCESS);
721 DNDS_assert(lapl_id >= 0);
722 if (coll_on_meta)
723 herr = H5Pset_all_coll_metadata_ops(lapl_id, offset.isDist() || offset == ArrayGlobalOffset_One), H5CHECK_Set;
724 htri_t exists_single = H5Lexists(group_id, name.c_str(), lapl_id);
725 htri_t exists_rank_offsets = H5Lexists(group_id, (name + "::rank_offsets").c_str(), lapl_id);
726 herr = H5Pclose(lapl_id), H5CHECK_Close;
727
728 if (exists_single > 0)
729 {
730 if (exists_rank_offsets > 0)
731 {
732 index nGlobal_offsets{-1};
733 int dim2_offsets{-1};
734 H5_ReadDataset(group_id, (name + "::rank_offsets").c_str(), nGlobal_offsets, -1, -1, DNDS_H5T_INDEX(), dxpl_id, dapl_id, nullptr, dim2_offsets);
735 DNDS_assert(dim2_offsets == -1);
736 if (nGlobal_offsets == mpi.size + 1)
737 {
738 std::array<index, 2> offsets = {-1, -1};
739 H5_ReadDataset(group_id, (name + "::rank_offsets").c_str(), nGlobal_offsets, mpi.rank, 2, DNDS_H5T_INDEX(), dxpl_id, dapl_id, offsets.data(), dim2_offsets);
740 index nGlobal{-1};
741 int dim2{-1};
742 H5_ReadDataset(group_id, name.c_str(), nGlobal, -1, -1, T_H5TYPE, dxpl_id, dapl_id, nullptr, dim2);
743 DNDS_assert(dim2 == -1);
744 DNDS_assert(offsets[1] >= offsets[0]);
745 offset = ArrayGlobalOffset{nGlobal, offsets[0]};
746 size = offsets[1] - offsets[0];
747 }
748 else
749 {
750 DNDS_assert_info(false, "no valid determination" +
751 fmt::format(" [{} {}], ", exists_single, exists_rank_offsets) + name);
752 }
753 }
754 else
755 {
756 DNDS_assert_info(false, "no valid determination; arrays of ArrayGlobalOffset_One need to be explicitly designated" +
757 fmt::format(" [{} {}], ", exists_single, exists_rank_offsets) + name);
758 }
759 }
760 else
761 DNDS_assert_info(false, "no valid determination; no dataset found " +
762 fmt::format(" [{} {}], ", exists_single, exists_rank_offsets) + name);
763 }
764 else
765 {
766 if (offset == ArrayGlobalOffset_One)
767 {
768 index nGlobal{-1};
769 int dim2{-1};
770 if (v != nullptr)
771 {
772 H5_ReadDataset(group_id, name.c_str(), nGlobal, 0, size, T_H5TYPE, dxpl_id, dapl_id, v, dim2);
773 DNDS_assert(dim2 == -1);
774 DNDS_assert(size == nGlobal);
775 }
776 else
777 {
778 H5_ReadDataset(group_id, name.c_str(), nGlobal, -1, -1, T_H5TYPE, dxpl_id, dapl_id, nullptr, dim2);
779 DNDS_assert(dim2 == -1);
780 size = nGlobal;
781 }
782 }
783 // else if (offset == ArrayGlobalOffset_Parts)
784 // {
785 // index nGlobal{-1};
786 // int dim2{-1};
787 // if (v != nullptr)
788 // {
789 // H5_ReadDataset(group_id, (name + "::" + std::to_string(mpi.rank)).c_str(), nGlobal, 0, size, T_H5TYPE, dxpl_id, dapl_id, v, dim2);
790 // DNDS_assert(dim2 == -1);
791 // DNDS_assert_info(size == nGlobal, fmt::format(", {}, {}", size, nGlobal));
792 // }
793 // else
794 // {
795 // H5_ReadDataset(group_id, (name + "::" + std::to_string(mpi.rank)).c_str(), nGlobal, -1, -1, T_H5TYPE, dxpl_id, dapl_id, nullptr, dim2);
796 // DNDS_assert(dim2 == -1);
797 // size = nGlobal;
798 // }
799 // }
800 else if (offset == ArrayGlobalOffset_EvenSplit)
801 {
802 // Even-split read: each rank reads ~N_global/nRanks rows.
803 // First pass (v==nullptr): query global size, compute local range, set offset.
804 // Second pass (v!=nullptr): read the data slice.
805 index nGlobal{-1};
806 int dim2{-1};
807 if (v == nullptr)
808 {
809 H5_ReadDataset(group_id, name.c_str(), nGlobal, -1, -1, T_H5TYPE, dxpl_id, dapl_id, nullptr, dim2);
810 DNDS_assert(dim2 == -1);
811 auto [start, end] = EvenSplitRange(mpi.rank, mpi.size, nGlobal);
812 size = end - start;
813 offset = ArrayGlobalOffset{index(size), start};
814 }
815 else
816 {
817 DNDS_assert_info(offset.isDist(), "EvenSplit second pass must have a resolved offset");
818 H5_ReadDataset(group_id, name.c_str(), nGlobal, offset.offset(), size, T_H5TYPE, dxpl_id, dapl_id, v, dim2);
819 DNDS_assert(dim2 == -1);
820 }
821 }
822 else if (offset.isDist())
823 {
824 if (v == nullptr)
825 {
826 // First pass (size query): local size is already in offset.size().
827 size = offset.size();
828 }
829 else
830 {
831 // Second pass (data read): read `size` elements at offset.offset().
832 // size==0 is valid for ranks with empty partitions (nGlobal < nRanks).
833 index nGlobal{-1};
834 int dim2{-1};
835 H5_ReadDataset(group_id, name.c_str(), nGlobal, offset.offset(), size, T_H5TYPE, dxpl_id, dapl_id, v, dim2);
836 DNDS_assert(dim2 == -1);
837 }
838 }
839 else
840 DNDS_assert_info(false, "offset ill-formed");
841 }
842
843 herr = H5Pclose(dxpl_id), H5CHECK_Close;
844 herr = H5Pclose(dapl_id), H5CHECK_Close;
845 herr = H5Gclose(group_id), H5CHECK_Close;
846 }
847
848 void SerializerH5::ReadIndexVector(const std::string &name, std::vector<index> &v, ArrayGlobalOffset &offset)
849 {
850 size_t size{0};
851 index dummy{};
852 ReadDataVector<index>(name, nullptr, size, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
853 v.resize(size);
854 DNDS_assert(!(offset == ArrayGlobalOffset_Unknown));
855 ReadDataVector<index>(name, size == 0 ? &dummy : v.data(), size, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
856 }
857 void SerializerH5::ReadRowsizeVector(const std::string &name, std::vector<rowsize> &v, ArrayGlobalOffset &offset)
858 {
859 size_t size{0};
860 rowsize dummy{};
861 ReadDataVector<rowsize>(name, nullptr, size, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
862 v.resize(size);
863 DNDS_assert(!(offset == ArrayGlobalOffset_Unknown));
864 ReadDataVector<rowsize>(name, size == 0 ? &dummy : v.data(), size, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
865 }
866 void SerializerH5::ReadRealVector(const std::string &name, std::vector<real> &v, ArrayGlobalOffset &offset)
867 {
868 size_t size{0};
869 real dummy{};
870 ReadDataVector<real>(name, nullptr, size, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
871 v.resize(size);
872 DNDS_assert(!(offset == ArrayGlobalOffset_Unknown));
873 ReadDataVector<real>(name, size == 0 ? &dummy : v.data(), size, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
874 }
875 void SerializerH5::ReadString(const std::string &name, std::string &v)
876 {
877 ReadAttributeScalar<std::string>(name, v, h5file, reading, cP, collectiveMetadataRW);
878 }
880 {
881 using tValue = host_device_vector<index>;
882 herr_t herr = 0;
883 std::string refPath;
884 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, collectiveMetadataRW);
885 htri_t exists_ref = H5Aexists(group_id, (name + "::ref").c_str());
886 herr = H5Gclose(group_id), H5CHECK_Close;
887 if (exists_ref > 0)
888 {
889 this->ReadString(name + "::ref", refPath);
890 }
891 else
892 {
893 refPath = cP + "/" + name;
894 }
895
896 if (pth_2_ssp.count(refPath))
897 {
898 // Dedup registry stores type-erased `ssp<tValue> *`; caller
899 // guarantees the stored type matches tValue.
900 v = *reinterpret_cast<ssp<tValue> *>(pth_2_ssp[refPath]);
901 }
902 else
903 {
904 v = std::make_shared<tValue>();
905 pth_2_ssp[refPath] = &v;
906
907 size_t size = 0;
908 index dummy{};
909 ReadDataVector<index>(refPath, nullptr, size, offset, h5file, reading, "/", mpi, collectiveMetadataRW, collectiveDataRW);
910 v->resize(size);
911 DNDS_assert(!(offset == ArrayGlobalOffset_Unknown));
912 ReadDataVector<index>(refPath, size == 0 ? &dummy : v->data(), size, offset, h5file, reading, "/", mpi, collectiveMetadataRW, collectiveDataRW);
913 }
914 }
916 {
917 using tValue = host_device_vector<rowsize>;
918 herr_t herr = 0;
919 std::string refPath;
920 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, collectiveMetadataRW);
921 htri_t exists_ref = H5Aexists(group_id, (name + "::ref").c_str());
922 herr = H5Gclose(group_id), H5CHECK_Close;
923 if (exists_ref > 0)
924 {
925 this->ReadString(name + "::ref", refPath);
926 }
927 else
928 {
929 refPath = cP + "/" + name;
930 }
931
932 if (pth_2_ssp.count(refPath))
933 {
934 // Dedup registry stores type-erased `ssp<tValue> *`; caller
935 // guarantees the stored type matches tValue.
936 v = *reinterpret_cast<ssp<tValue> *>(pth_2_ssp[refPath]);
937 }
938 else
939 {
940 v = std::make_shared<tValue>();
941 pth_2_ssp[refPath] = &v;
942
943 size_t size = 0;
944 rowsize dummy{};
945 ReadDataVector<rowsize>(refPath, nullptr, size, offset, h5file, reading, "/", mpi, collectiveMetadataRW, collectiveDataRW);
946 v->resize(size);
947 DNDS_assert(!(offset == ArrayGlobalOffset_Unknown));
948 ReadDataVector<rowsize>(refPath, size == 0 ? &dummy : v->data(), size, offset, h5file, reading, "/", mpi, collectiveMetadataRW, collectiveDataRW);
949 }
950 }
951 void SerializerH5::WriteUint8Array(const std::string &name, const uint8_t *data, index size, ArrayGlobalOffset offset)
952 {
953 WriteDataVector<uint8_t>(name, data, size, offset, chunksize, deflateLevel, h5file, reading, cP, mpi, commDup, collectiveMetadataRW, collectiveDataRW);
954 }
955 void SerializerH5::ReadUint8Array(const std::string &name, uint8_t *data, index &size, ArrayGlobalOffset &offset)
956 {
957 size_t size_size_t{0};
958 if (data == nullptr)
959 {
960 ReadDataVector<uint8_t>(name, nullptr, size_size_t, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
961 size = size_size_t; // todo: check overflow?
962 DNDS_assert(!(offset == ArrayGlobalOffset_Unknown));
963 }
964 else
965 {
966 DNDS_assert(!(offset == ArrayGlobalOffset_Unknown));
967 size_size_t = size;
968 ReadDataVector<uint8_t>(name, data, size_size_t, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
969 }
970 }
971}
#define DNDS_assert_info(expr, info)
Debug-only assertion with an extra std::string info message.
Definition Errors.hpp:117
#define DNDS_assert(expr)
Debug-only assertion (compiled out when DNDS_NDEBUG is defined). Prints the expression + file/line + ...
Definition Errors.hpp:112
HDF5 native type mappings for DNDS index, rowsize, and real types.
std::vector< int64_t > offsets
Per-cell cumulative offset into connectivity.
#define H5CHECK_Iter
#define H5CHECK_Set
#define H5CHECK_Close
MPI-parallel HDF5 serializer implementing the SerializerBase interface.
Describes one rank's window into a globally-distributed dataset.
void dedupRegister(const ssp< T > &v, const std::string &path)
Register a shared pointer after writing its data.
bool dedupLookup(const ssp< T > &v, std::string &outPath)
Check if a shared pointer was already written; if so return its path.
std::map< std::string, void * > pth_2_ssp
Reverse map for read-side dedup: path -> raw pointer to the ssp local variable.
void dedupClear()
Clear all dedup state (call on CloseFile).
void WriteIndex(const std::string &name, index v) override
Write a scalar index under name.
void ReadRealVector(const std::string &name, std::vector< real > &v, ArrayGlobalOffset &offset) override
std::string GetCurrentPath() override
String form of the current path.
void WriteSharedIndexVector(const std::string &name, const ssp< host_device_vector< index > > &v, ArrayGlobalOffset offset) override
Write a shared index vector; deduplicated across multiple writes that share the same shared_ptr.
void ReadSharedRowsizeVector(const std::string &name, ssp< host_device_vector< rowsize > > &v, ArrayGlobalOffset &offset) override
void ReadIndexVector(const std::string &name, std::vector< index > &v, ArrayGlobalOffset &offset) override
void ReadIndex(const std::string &name, index &v) override
Read a scalar index into v.
void WriteInt(const std::string &name, int v) override
Write a scalar int under name at the current path.
void WriteUint8Array(const std::string &name, const uint8_t *data, index size, ArrayGlobalOffset offset) override
Write a raw byte buffer under name. offset.isDist() = true means the caller provides the exact per-ra...
std::set< std::string > ListCurrentPath() override
Names of direct children of the current path.
void ReadInt(const std::string &name, int &v) override
Read a scalar int into v.
void OpenFile(const std::string &fName, bool read) override
Open a backing file (H5 file or JSON file depending on subclass).
void WriteRealVector(const std::string &name, const std::vector< real > &v, ArrayGlobalOffset offset) override
Write a real vector (collective for H5).
void GoToPath(const std::string &p) override
Navigate to an existing path. Supports / -separated segments.
void WriteReal(const std::string &name, real v) override
Write a scalar real under name.
void ReadReal(const std::string &name, real &v) override
Read a scalar real into v.
void ReadRowsizeVector(const std::string &name, std::vector< rowsize > &v, ArrayGlobalOffset &offset) override
void CloseFile() override
Close the backing file, flushing buffers.
void WriteString(const std::string &name, const std::string &v) override
Write a UTF-8 string under name.
void CreatePath(const std::string &p) override
Create a sub-path (H5 group / JSON object) at the current location.
void ReadString(const std::string &name, std::string &v) override
Read a UTF-8 string into v.
void WriteSharedRowsizeVector(const std::string &name, const ssp< host_device_vector< rowsize > > &v, ArrayGlobalOffset offset) override
Write a shared rowsize vector; deduplicated across multiple writes.
void WriteRowsizeVector(const std::string &name, const std::vector< rowsize > &v, ArrayGlobalOffset offset) override
Write a rowsize vector (collective for H5).
void ReadSharedIndexVector(const std::string &name, ssp< host_device_vector< index > > &v, ArrayGlobalOffset &offset) override
void WriteIndexVector(const std::string &name, const std::vector< index > &v, ArrayGlobalOffset offset) override
Write an index vector (collective for H5). offset carries the distribution mode (ArrayGlobalOffset_Pa...
void ReadUint8Array(const std::string &name, uint8_t *data, index &size, ArrayGlobalOffset &offset) override
Two-pass byte array read.
Host + optional device vector of trivially copyable T.
Definition Vector.hpp:217
MPI_int Bcast(void *buf, MPI_int num, MPI_Datatype type, MPI_int source_rank, MPI_Comm comm)
dumb wrapper
Definition MPI.cpp:149
MPI_int Scan(const void *sendbuf, void *recvbuf, MPI_int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
Wrapper over MPI_Scan (inclusive prefix reduction).
Definition MPI.cpp:219
void AllreduceOneIndex(index &v, MPI_Op op, const MPIInfo &mpi)
Single-scalar Allreduce helper for indices (in-place, count = 1).
Definition MPI.hpp:727
void WriteAttributeScalar(const std::string &name, const T &v, hid_t h5file, bool reading, const std::string &cP, bool coll_on_meta)
void ReadAttributeScalar(const std::string &name, T &v, hid_t h5file, bool reading, const std::string &cP, bool coll_on_meta)
hid_t DNDS_H5T_REAL()
HDF5 native datatype matching DNDS real (currently H5T_NATIVE_DOUBLE).
Definition HDF5.hpp:36
int32_t rowsize
Row-width / per-row element-count type (signed 32-bit).
Definition Defines.hpp:114
std::pair< index, index > EvenSplitRange(int rank, int nRanks, index nGlobal)
Split a global range [0, nGlobal) evenly among nRanks workers.
Definition Defines.hpp:134
std::vector< std::string > splitSString(const std::string &str, char delim)
Definition Defines.hpp:777
int64_t index
Global row / DOF index type (signed 64-bit; handles multi-billion-cell meshes).
Definition Defines.hpp:112
std::shared_ptr< T > ssp
Shortened alias for std::shared_ptr used pervasively in DNDSR.
Definition Defines.hpp:143
hid_t DNDS_H5T_ROWSIZE()
HDF5 native datatype matching DNDS rowsize (currently H5T_NATIVE_INT32).
Definition HDF5.hpp:24
double real
Canonical floating-point scalar used throughout DNDSR (double precision).
Definition Defines.hpp:110
hid_t DNDS_H5T_INDEX()
HDF5 native datatype matching DNDS index (currently H5T_NATIVE_INT64).
Definition HDF5.hpp:12
std::ostream & log()
Return the current DNDSR log stream (either std::cout or the installed file).
Definition Defines.cpp:50
Lightweight bundle of an MPI communicator and the calling rank's coordinates.
Definition MPI.hpp:231
std::vector< std::string > attributes
std::vector< std::string > groups
std::vector< std::string > datasets
Eigen::Matrix< real, 5, 1 > v
constexpr DNDS::index nLocal
const tPoint const tPoint const tPoint & p