DNDSR 0.1.0.dev1+gcd065ad
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>
5#include "DNDS/HDF5.hpp"
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 {
181 std::string current_path;
183 std::string get_indent() const
184 {
185 return std::string(std::count(current_path.begin(), current_path.end(), '/') * 2, ' ');
186 }
187 };
188
189 static herr_t link_iterate_cb(hid_t group_id, const char *name, const H5L_info_t *info, void *op_data)
190 {
191 herr_t herr;
192 auto *data = static_cast<TraverseData *>(op_data);
193 bool coll_on_meta = data->coll_on_meta;
194 std::string full_name = data->current_path + "/" + name;
195 // We only care about hard links that represent HDF5 objects
196 // std::cout << "name is " << name << std::endl;
197 if (info->type == H5L_TYPE_HARD)
198 {
199 data->contents.groups.push_back(name);
200 return 0;
201 //! now it seems obj_info is not correctly retrieved with type=Unknown
202 //! TODO: fix this and get actual object type
203 H5O_info_t obj_info;
204 hid_t lapl_id = H5Pcreate(H5P_LINK_ACCESS);
205 DNDS_assert(lapl_id >= 0);
206 herr = H5Pset_all_coll_metadata_ops(lapl_id, true), H5CHECK_Set;
207 // if (coll_on_meta)
208 // herr = H5Pset_all_coll_metadata_ops(lapl_id, true), H5CHECK_Set;
209 if (H5Oget_info_by_name(group_id, name, &obj_info, H5P_DEFAULT, lapl_id) >= 0)
210 {
211 log() << "name is " << name << ", \n"
212 << obj_info.type << "\n"
213 << obj_info.mtime << "\n"
214 << obj_info.rc << "\n"
215 << obj_info.fileno << "\n"
216 << std::endl;
217 switch (obj_info.type)
218 {
219 case H5O_TYPE_GROUP:
220 {
221 data->contents.groups.push_back(name);
222 // std::cout << data->get_indent() << "Group: " << full_name << std::endl;
223 break;
224 }
225 case H5O_TYPE_DATASET:
226 {
227 data->contents.datasets.push_back(name);
228 // std::cout << data->get_indent() << "Dataset: " << full_name << std::endl;
229 break;
230 }
231 default:
232 // Other types (e.g., symbols, but less common for direct user interaction)
233 break;
234 }
235 }
236 else
237 {
238 log() << "Error getting object info for: " << full_name << std::endl;
239 H5Eprint2(H5E_DEFAULT, stderr);
240 }
241 herr = H5Pclose(lapl_id), H5CHECK_Close;
242 }
243 else
244 {
245 // Handle soft links or external links if their names are also desired
246 // std::cout << data->get_indent() << "Link: " << full_name << " (Non-hard link)" << std::endl;
247 }
248 return 0; // Continue iteration
249 }
250
251 // Callback for H5Aiterate (iterating attributes)
252 static herr_t attribute_iterate_cb(hid_t obj_id, const char *attr_name, const H5A_info_t *info, void *op_data)
253 {
254 TraverseData *data = static_cast<TraverseData *>(op_data);
255 data->contents.attributes.push_back(attr_name);
256 // std::string full_attr_name = data->current_path + "@" + attr_name; // Common convention for attribute paths
257 // std::cout << data->get_indent() << " Attribute: " << full_attr_name << std::endl; // Indent more for attributes
258 return 0; // Continue iteration
259 }
260
261 std::set<std::string> SerializerH5::ListCurrentPath()
262 {
263 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, collectiveMetadataRW);
264 H5Contents contents;
265 TraverseData data{contents, cP, collectiveMetadataRW};
266 herr_t herr{0};
267 herr = H5Aiterate(group_id, H5_INDEX_NAME, H5_ITER_INC, NULL, attribute_iterate_cb, &data), H5CHECK_Iter;
268 herr = H5Literate(group_id, H5_INDEX_NAME, H5_ITER_INC, NULL, link_iterate_cb, &data), H5CHECK_Iter;
269 H5Gclose(group_id), H5CHECK_Close; // don't forget this
270 std::set<std::string> ret;
271 for (auto &v : contents.attributes)
272 ret.insert(v);
273 for (auto &v : contents.groups)
274 ret.insert(v);
275 for (auto &v : contents.datasets)
276 ret.insert(v);
277 return ret;
278 }
279
280 template <typename T>
281 // using T = int;
282 void WriteAttributeScalar(const std::string &name, const T &v,
283 hid_t h5file, bool reading, const std::string &cP,
284 bool coll_on_meta)
285 {
286 hid_t T_H5TYPE = -1;
287 if constexpr (std::is_same_v<T, int>)
288 T_H5TYPE = H5T_NATIVE_INT;
289 else if constexpr (std::is_same_v<T, index>)
290 T_H5TYPE = DNDS_H5T_INDEX();
291 else if constexpr (std::is_same_v<T, rowsize>)
292 T_H5TYPE = DNDS_H5T_ROWSIZE();
293 else if constexpr (std::is_same_v<T, real>)
294 T_H5TYPE = DNDS_H5T_REAL();
295 else if constexpr (std::is_same_v<T, std::string>)
296 ; // pass
297 else
298 static_assert(std::is_same_v<T, real>);
299
300 T vV = v;
301
302 if constexpr (!std::is_same_v<T, std::string>)
303 {
304 herr_t herr{0};
305 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, coll_on_meta);
306 hsize_t attr_size = 1;
307 hid_t attr_space = H5Screate(H5S_SCALAR);
308 hid_t aapl_id = H5Pcreate(H5P_ATTRIBUTE_ACCESS);
309 if (coll_on_meta)
310 herr = H5Pset_all_coll_metadata_ops(aapl_id, true), H5CHECK_Set;
311 hid_t attr_id = H5Acreate(group_id, name.c_str(), T_H5TYPE, attr_space, H5P_DEFAULT, aapl_id);
312 herr = H5Awrite(attr_id, T_H5TYPE, &vV), H5CHECK_Set;
313 H5Aclose(attr_id);
314 H5Pclose(aapl_id);
315 H5Sclose(attr_space);
316 H5Gclose(group_id);
317 }
318 else
319 {
320 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, coll_on_meta);
321
322 herr_t herr{0};
323 hid_t aapl_id = H5Pcreate(H5P_ATTRIBUTE_ACCESS);
324 if (coll_on_meta)
325 herr = H5Pset_all_coll_metadata_ops(aapl_id, true), H5CHECK_Set;
326 hid_t scalar_space = H5Screate(H5S_SCALAR);
327 hid_t string_type = H5Tcreate(H5T_STRING, v.length());
328 herr = H5Tset_strpad(string_type, H5T_STR_NULLPAD), H5CHECK_Set;
329 hid_t type_attr_id = H5Acreate(group_id, name.c_str(), string_type, scalar_space, H5P_DEFAULT, aapl_id);
330 DNDS_assert_info(type_attr_id >= 0, fmt::format("{}, {}, {}", name, v, cP));
331 herr = H5Awrite(type_attr_id, string_type, v.data()), H5CHECK_Set;
332
333 herr = H5Aclose(type_attr_id), H5CHECK_Close;
334 herr = H5Tclose(string_type), H5CHECK_Close;
335 herr = H5Pclose(aapl_id), H5CHECK_Close;
336 herr = H5Sclose(scalar_space), H5CHECK_Close;
337 herr = H5Gclose(group_id), H5CHECK_Close;
338 }
339 }
340
341 void SerializerH5::WriteInt(const std::string &name, int v)
342 {
343 WriteAttributeScalar<int>(name, v, h5file, reading, cP, collectiveMetadataRW);
344 }
345 void SerializerH5::WriteIndex(const std::string &name, index v)
346 {
347 WriteAttributeScalar<index>(name, v, h5file, reading, cP, collectiveMetadataRW);
348 }
349 void SerializerH5::WriteReal(const std::string &name, real v)
350 {
351 WriteAttributeScalar<real>(name, v, h5file, reading, cP, collectiveMetadataRW);
352 }
353
354 static void H5_WriteDataset(hid_t loc, const char *name, index nGlobal, index nOffset, index nLocal,
355 hid_t file_dataType, hid_t mem_dataType, hid_t dxpl_id, hid_t dapl_id, int64_t chunksize, int deflateLevel,
356 const void *buf, int dim2 = -1)
357 {
358 int herr{0};
359 DNDS_assert_info(nGlobal >= 0 && nLocal >= 0 && nOffset >= 0,
360 fmt::format("{},{},{}", nGlobal, nLocal, nOffset));
361 int rank = dim2 >= 0 ? 2 : 1;
362 std::array<hsize_t, 2> ranksFull{hsize_t(nGlobal), hsize_t(dim2)};
363 std::array<hsize_t, 2> ranksFullUnlim{chunksize > 0 ? H5S_UNLIMITED : hsize_t(nGlobal), hsize_t(dim2)};
364 std::array<hsize_t, 2> offset{hsize_t(nOffset), 0};
365 std::array<hsize_t, 2> siz{hsize_t(nLocal), hsize_t(dim2)};
366 hid_t memSpace = H5Screate_simple(rank, siz.data(), NULL);
367 hid_t fileSpace = H5Screate_simple(rank, ranksFull.data(), ranksFullUnlim.data());
368 std::array<hsize_t, 2> chunk_dims{hsize_t(chunksize > 0 ? chunksize : 0), dim2 >= 0 ? hsize_t(dim2) : 0};
369 hid_t dcpl_id = H5Pcreate(H5P_DATASET_CREATE);
370 if (chunk_dims[0] > 0)
371 herr = H5Pset_chunk(dcpl_id, dim2 > 0 ? 2 : 1, chunk_dims.data()), H5CHECK_Set;
372#ifdef H5_HAVE_FILTER_DEFLATE
373 if (deflateLevel > 0)
374 herr = H5Pset_deflate(dcpl_id, deflateLevel), H5CHECK_Set;
375#endif
376
377 hid_t dset_id = H5Dcreate(loc, name, file_dataType, fileSpace, H5P_DEFAULT, dcpl_id, dapl_id);
378 DNDS_assert_info(H5I_INVALID_HID != dset_id, "dataset create failed");
379 herr = H5Sclose(fileSpace);
380 fileSpace = H5Dget_space(dset_id);
381 herr |= H5Sselect_hyperslab(fileSpace, H5S_SELECT_SET, offset.data(), NULL, siz.data(), NULL);
382 herr |= H5Dwrite(dset_id, mem_dataType, memSpace, fileSpace, dxpl_id, buf);
383 herr |= H5Dclose(dset_id);
384 herr |= H5Pclose(dcpl_id);
385 herr |= H5Sclose(fileSpace);
386 herr |= H5Sclose(memSpace);
387 DNDS_assert_info(herr >= 0,
388 "h5 error " + fmt::format(
389 "nGlobal {}, nOffset {}, nLocal {}, name {}",
390 nGlobal, nOffset, nLocal, name));
391 } //! currently a direct copy
392
393 template <typename T = index>
394 // using T = index;
395 static void WriteDataVector(const std::string &name, const T *v, size_t size, ArrayGlobalOffset offset, int64_t chunksize, int deflateLevel,
396 hid_t h5file, bool reading, const std::string &cP, const MPIInfo &mpi, MPI_Comm commDup,
397 bool coll_on_meta, bool coll_on_data)
398 {
399 hid_t T_H5TYPE = H5T_NATIVE_INT;
400 if constexpr (std::is_same_v<T, index>)
401 T_H5TYPE = DNDS_H5T_INDEX();
402 else if constexpr (std::is_same_v<T, rowsize>)
403 T_H5TYPE = DNDS_H5T_ROWSIZE();
404 else if constexpr (std::is_same_v<T, real>)
405 T_H5TYPE = DNDS_H5T_REAL();
406 else if constexpr (std::is_same_v<T, uint8_t>)
407 T_H5TYPE = H5T_NATIVE_UINT8;
408 else
409 static_assert(std::is_same_v<T, uint8_t>);
410
411 herr_t herr{0};
412 hid_t dxpl_id = H5Pcreate(H5P_DATASET_XFER);
413 //! is this necessary? seems even global-unique vector needs this if deflate is on!
414 if (((offset.isDist() || offset == ArrayGlobalOffset_Parts) && (coll_on_data)) || deflateLevel > 0)
415 herr = H5Pset_dxpl_mpio(dxpl_id, H5FD_MPIO_COLLECTIVE), H5CHECK_Set;
416 else
417 herr = H5Pset_dxpl_mpio(dxpl_id, H5FD_MPIO_INDEPENDENT), H5CHECK_Set;
418 DNDS_assert_info(herr >= 0, "h5 error");
419 // if is dist array, we use coll metadata
420 hid_t dapl_id = H5Pcreate(H5P_DATASET_ACCESS);
421 if (coll_on_meta)
422 herr = H5Pset_all_coll_metadata_ops(dapl_id, true), H5CHECK_Set;
423
424 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, coll_on_meta);
425 if (offset == ArrayGlobalOffset_One)
426 H5_WriteDataset(group_id, name.c_str(), size, 0, mpi.rank == 0 ? size : 0,
427 T_H5TYPE, T_H5TYPE, dxpl_id, dapl_id, chunksize, deflateLevel,
428 v);
429 else if (offset == ArrayGlobalOffset_Parts) // now we force it to be
430 {
431 uint64_t sizeU{size}, sizeOff{0}, sizeGlobal{0};
432 MPI::Scan(&sizeU, &sizeOff, 1, MPI_UINT64_T, MPI_SUM, commDup);
433 sizeGlobal = sizeOff;
434 MPI::Bcast(&sizeGlobal, 1, MPI_UINT64_T, mpi.size - 1, commDup);
435 sizeOff -= sizeU;
436
437 H5_WriteDataset(group_id, name.c_str(), sizeGlobal, sizeOff, size,
438 T_H5TYPE, T_H5TYPE, dxpl_id, dapl_id, chunksize, deflateLevel,
439 v);
440 // rank offset array
441 std::array<index, 2> offsetC = {index(sizeOff), index(sizeGlobal)};
442 H5_WriteDataset(group_id, (name + "::rank_offsets").c_str(), mpi.size + 1, mpi.rank, (mpi.rank == mpi.size - 1) ? 2 : 1,
443 DNDS_H5T_INDEX(), DNDS_H5T_INDEX(), dxpl_id, dapl_id, chunksize, deflateLevel,
444 offsetC.data());
445 }
446 else if (offset.isDist())
447 {
448
449 H5_WriteDataset(group_id, name.c_str(), offset.size(), offset.offset(), size,
450 T_H5TYPE, T_H5TYPE, dxpl_id, dapl_id, chunksize, deflateLevel,
451 v);
452 // rank offset array
453 std::array<index, 2> offsetC = {offset.offset(), offset.size()};
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
459 DNDS_assert_info(false, "offset ill-formed");
460
461 herr = H5Pclose(dxpl_id), H5CHECK_Close;
462 herr = H5Pclose(dapl_id), H5CHECK_Close;
463 herr = H5Gclose(group_id), H5CHECK_Close;
464 }
465
466 void SerializerH5::WriteIndexVector(const std::string &name, const std::vector<index> &v, ArrayGlobalOffset offset)
467 {
468 WriteDataVector<index>(name, v.data(), v.size(), offset, chunksize, deflateLevel, h5file, reading, cP, mpi, commDup, collectiveMetadataRW, collectiveDataRW);
469 }
470 void SerializerH5::WriteRowsizeVector(const std::string &name, const std::vector<rowsize> &v, ArrayGlobalOffset offset)
471 {
472 WriteDataVector<rowsize>(name, v.data(), v.size(), offset, chunksize, deflateLevel, h5file, reading, cP, mpi, commDup, collectiveMetadataRW, collectiveDataRW);
473 }
474 void SerializerH5::WriteRealVector(const std::string &name, const std::vector<real> &v, ArrayGlobalOffset offset)
475 {
476 WriteDataVector<real>(name, v.data(), v.size(), offset, chunksize, deflateLevel, h5file, reading, cP, mpi, commDup, collectiveMetadataRW, collectiveDataRW);
477 }
478 void SerializerH5::WriteString(const std::string &name, const std::string &v)
479 {
480 WriteAttributeScalar<std::string>(name, v, h5file, reading, cP, collectiveMetadataRW);
481 }
483 {
484 std::string refPath;
485 index CanNotShare = dedupLookup(v, refPath) ? 0 : 1;
486 MPI::AllreduceOneIndex(CanNotShare, MPI_MAX, MPIInfo{commDup, mpi.rank, mpi.size});
487 if (!CanNotShare)
488 this->WriteString(name + "::ref", refPath);
489 else
490 {
491 WriteDataVector<index>(name, v->data(), v->size(), offset, chunksize, deflateLevel, h5file, reading, cP, mpi, commDup, collectiveMetadataRW, collectiveDataRW);
492 dedupRegister(v, cP + "/" + name);
493 }
494 }
496 {
497 std::string refPath;
498 index CanNotShare = dedupLookup(v, refPath) ? 0 : 1;
499 MPI::AllreduceOneIndex(CanNotShare, MPI_MAX, MPIInfo{commDup, mpi.rank, mpi.size});
500 if (!CanNotShare)
501 this->WriteString(name + "::ref", refPath);
502 else
503 {
504 WriteDataVector<rowsize>(name, v->data(), v->size(), offset, chunksize, deflateLevel, h5file, reading, cP, mpi, commDup, collectiveMetadataRW, collectiveDataRW);
505 dedupRegister(v, cP + "/" + name);
506 }
507 }
508
509 template <typename T>
510 // using T = int;
511 void ReadAttributeScalar(const std::string &name, T &v,
512 hid_t h5file, bool reading, const std::string &cP,
513 bool coll_on_meta)
514 {
515 hid_t T_H5TYPE = -1;
516 if constexpr (std::is_same_v<T, int>)
517 T_H5TYPE = H5T_NATIVE_INT;
518 else if constexpr (std::is_same_v<T, index>)
519 T_H5TYPE = DNDS_H5T_INDEX();
520 else if constexpr (std::is_same_v<T, rowsize>)
521 T_H5TYPE = DNDS_H5T_ROWSIZE();
522 else if constexpr (std::is_same_v<T, real>)
523 T_H5TYPE = DNDS_H5T_REAL();
524 else if constexpr (std::is_same_v<T, std::string>)
525 ; // pass
526 else
527 static_assert(std::is_same_v<T, real>);
528
529 T vV = v;
530
531 if constexpr (!std::is_same_v<T, std::string>)
532 {
533 herr_t herr{0};
534 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, coll_on_meta);
535 hid_t aapl_id = H5Pcreate(H5P_ATTRIBUTE_ACCESS);
536 herr = H5Pset_all_coll_metadata_ops(aapl_id, true), H5CHECK_Set;
537 // herr = H5Pset_coll_metadata_write(aapl_id, true), H5CHECK_Set; // setting to present in dapl
538 hid_t attr_id = H5Aopen(group_id, name.c_str(), aapl_id);
539 DNDS_assert_info(attr_id >= 0, fmt::format("attempting to open attribute [{}/{}] failed", cP, name));
540 hid_t attr_space = H5Aget_space(attr_id);
541 int ndims = H5Sget_simple_extent_ndims(attr_space);
542 DNDS_assert_info(ndims == 0, fmt::format("attempting to read attribute [{}/{}] which is not scalar", cP, name));
543 herr = H5Aread(attr_id, T_H5TYPE, &v), H5CHECK_Set;
544 herr = H5Aclose(attr_id), H5CHECK_Close;
545 herr = H5Pclose(aapl_id), H5CHECK_Close;
546 herr = H5Sclose(attr_space), H5CHECK_Close;
547 herr = H5Gclose(group_id), H5CHECK_Close;
548 }
549 else
550 {
551 herr_t herr{0};
552 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, coll_on_meta);
553 hid_t aapl_id = H5Pcreate(H5P_ATTRIBUTE_ACCESS);
554 if (coll_on_meta)
555 {
556 herr = H5Pset_all_coll_metadata_ops(aapl_id, true), H5CHECK_Set;
557 // herr = H5Pset_coll_metadata_write(aapl_id, true), H5CHECK_Set; // setting not present in aapl
558 }
559 hid_t attr_id = H5Aopen(group_id, name.c_str(), aapl_id);
560 DNDS_assert_info(attr_id >= 0, fmt::format("attempting to open attribute [{}/{}] failed", cP, name));
561 hid_t attr_space = H5Aget_space(attr_id);
562 int ndims = H5Sget_simple_extent_ndims(attr_space);
563 DNDS_assert_info(ndims == 0, fmt::format("attempting to read attribute [{}/{}] which is not scalar", cP, name));
564 hid_t dtype_id = H5Aget_type(attr_id);
565 bool is_varlen = H5Tis_variable_str(dtype_id);
566
567 if (is_varlen)
568 {
569 // Variable-length string: HDF5 will allocate memory
570 char *attr_value = nullptr;
571 herr = H5Aread(attr_id, dtype_id, &attr_value), H5CHECK_Set;
572 // std::cout << "Read Attribute (Variable-Length): " << attr_value << "\n";
573 v = attr_value; // copy as null-terminated string
574 H5free_memory(attr_value); // Must free memory manually
575 }
576 else
577 {
578 // Fixed-length string: Allocate buffer based on datatype size
579 size_t size = H5Tget_size(dtype_id);
580 std::vector<char> buffer(size + 1, '\0'); // +1 for null terminator
581 herr = H5Aread(attr_id, dtype_id, buffer.data()), H5CHECK_Set;
582 // std::cout << "Read Attribute (Fixed-Length): " << buffer.data() << "\n";
583 v = buffer.data(); // copy as null-terminated string
584 }
585
586 herr = H5Tclose(dtype_id), H5CHECK_Close;
587 herr = H5Aclose(attr_id), H5CHECK_Close;
588 herr = H5Pclose(aapl_id), H5CHECK_Close;
589 herr = H5Sclose(attr_space), H5CHECK_Close;
590 herr = H5Gclose(group_id), H5CHECK_Close;
591 }
592 }
593
594 void SerializerH5::ReadInt(const std::string &name, int &v)
595 {
596 ReadAttributeScalar<int>(name, v, h5file, reading, cP, collectiveMetadataRW);
597 }
598 void SerializerH5::ReadIndex(const std::string &name, index &v)
599 {
600 ReadAttributeScalar<index>(name, v, h5file, reading, cP, collectiveMetadataRW);
601 }
602 void SerializerH5::ReadReal(const std::string &name, real &v)
603 {
604 ReadAttributeScalar<real>(name, v, h5file, reading, cP, collectiveMetadataRW);
605 }
606
607 /**
608 * @brief
609 *
610 * when #buf == 2, only set nGlobal and dim2
611 *
612 * @param loc
613 * @param name
614 * @param nGlobal
615 * @param nOffset
616 * @param nLocal
617 * @param file_dataType
618 * @param mem_dataType
619 * @param plist_id
620 * @param dcpl_id
621 * @param dapl_id
622 * @param buf
623 * @param dim2
624 */
625 // When dxpl_id uses collective I/O (H5FD_MPIO_COLLECTIVE), every rank MUST
626 // call H5Dread even if nLocal==0. Callers must provide a valid (possibly
627 // dummy) buf pointer; do NOT skip this function based on nLocal==0.
628 static void H5_ReadDataset(hid_t loc, const char *name, index &nGlobal, index nOffset, index nLocal,
629 hid_t mem_dataType, hid_t dxpl_id, hid_t dapl_id,
630 void *buf, int &dim2)
631 {
632 int herr{0};
633 hid_t dset_id = H5Dopen(loc, name, dapl_id);
634 DNDS_assert_info(dset_id >= 0, fmt::format("dataset [{}] open failed", name));
635 hid_t fileSpace = H5Dget_space(dset_id);
636 DNDS_assert_info(fileSpace >= 0, fmt::format("dataset [{}] filespace open failed", name));
637 int ndims = H5Sget_simple_extent_ndims(fileSpace);
638 DNDS_assert_info(ndims == 1 || ndims == 2, fmt::format("dataset [{}] not having 1 or 2 dims!", name));
639 std::array<hsize_t, 2> sizes;
640 ndims = H5Sget_simple_extent_dims(fileSpace, sizes.data(), nullptr);
641 if (ndims == 2)
642 dim2 = sizes[1];
643 else
644 dim2 = -1;
645 nGlobal = sizes[0];
646
647 if (buf != nullptr)
648 {
649 int rank = dim2 >= 0 ? 2 : 1;
650 std::array<hsize_t, 2> offset{hsize_t(nOffset), 0};
651 std::array<hsize_t, 2> siz{hsize_t(nLocal), hsize_t(dim2)};
652 hid_t memSpace = H5Screate_simple(rank, siz.data(), NULL);
653 DNDS_assert(memSpace > 0);
654 herr = H5Sselect_hyperslab(fileSpace, H5S_SELECT_SET, offset.data(), NULL, siz.data(), NULL), H5CHECK_Set;
655 herr = H5Dread(dset_id, mem_dataType, memSpace, fileSpace, dxpl_id, buf), H5CHECK_Set;
656 herr = H5Sclose(memSpace), H5CHECK_Close;
657 }
658 herr = H5Dclose(dset_id), H5CHECK_Close;
659 herr = H5Sclose(fileSpace), H5CHECK_Close;
660 }
661
662 template <typename T = index>
663 // using T = index;
664 static void ReadDataVector(const std::string &name, T *v, size_t &size, ArrayGlobalOffset &offset,
665 hid_t h5file, bool reading, const std::string &cP, const MPIInfo &mpi,
666 bool coll_on_meta, bool coll_on_data)
667 {
668 hid_t T_H5TYPE = H5T_NATIVE_INT;
669 if constexpr (std::is_same_v<T, index>)
670 T_H5TYPE = DNDS_H5T_INDEX();
671 else if constexpr (std::is_same_v<T, rowsize>)
672 T_H5TYPE = DNDS_H5T_ROWSIZE();
673 else if constexpr (std::is_same_v<T, real>)
674 T_H5TYPE = DNDS_H5T_REAL();
675 else if constexpr (std::is_same_v<T, uint8_t>)
676 T_H5TYPE = H5T_NATIVE_UINT8;
677 else
678 static_assert(std::is_same_v<T, uint8_t>);
679
680 herr_t herr{0};
681 hid_t dxpl_id = H5Pcreate(H5P_DATASET_XFER);
682 if ((offset.isDist() || offset == ArrayGlobalOffset_One) && coll_on_data) //! is this necessary?
683 herr = H5Pset_dxpl_mpio(dxpl_id, H5FD_MPIO_COLLECTIVE), H5CHECK_Set;
684 else
685 herr = H5Pset_dxpl_mpio(dxpl_id, H5FD_MPIO_INDEPENDENT), H5CHECK_Set;
686
687 // if is dist array, we use coll metadata
688 hid_t dapl_id = H5Pcreate(H5P_DATASET_ACCESS);
689 if (coll_on_meta)
690 {
691 herr = H5Pset_all_coll_metadata_ops(dapl_id, offset.isDist() || offset == ArrayGlobalOffset_One), H5CHECK_Set;
692 // herr = H5Pset_coll_metadata_write(dapl_id, offset.isDist() || offset == ArrayGlobalOffset_One), H5CHECK_Set; // tsetting to present in dapl
693 }
694 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, coll_on_meta);
695
696 if (offset == ArrayGlobalOffset_Unknown) // need detection, auto ArrayGlobalOffset_Parts or distributed array
697 {
698 DNDS_assert(v == nullptr); // must be a size-query call
699
700 hid_t lapl_id = H5Pcreate(H5P_LINK_ACCESS);
701 DNDS_assert(lapl_id >= 0);
702 if (coll_on_meta)
703 herr = H5Pset_all_coll_metadata_ops(lapl_id, offset.isDist() || offset == ArrayGlobalOffset_One), H5CHECK_Set;
704 htri_t exists_single = H5Lexists(group_id, name.c_str(), lapl_id);
705 htri_t exists_rank_offsets = H5Lexists(group_id, (name + "::rank_offsets").c_str(), lapl_id);
706 herr = H5Pclose(lapl_id), H5CHECK_Close;
707
708 if (exists_single > 0)
709 {
710 if (exists_rank_offsets > 0)
711 {
712 index nGlobal_offsets{-1};
713 int dim2_offsets{-1};
714 H5_ReadDataset(group_id, (name + "::rank_offsets").c_str(), nGlobal_offsets, -1, -1, DNDS_H5T_INDEX(), dxpl_id, dapl_id, nullptr, dim2_offsets);
715 DNDS_assert(dim2_offsets == -1);
716 if (nGlobal_offsets == mpi.size + 1)
717 {
718 std::array<index, 2> offsets = {-1, -1};
719 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);
720 index nGlobal{-1};
721 int dim2{-1};
722 H5_ReadDataset(group_id, name.c_str(), nGlobal, -1, -1, T_H5TYPE, dxpl_id, dapl_id, nullptr, dim2);
723 DNDS_assert(dim2 == -1);
724 DNDS_assert(offsets[1] >= offsets[0]);
725 offset = ArrayGlobalOffset{nGlobal, offsets[0]};
726 size = offsets[1] - offsets[0];
727 }
728 else
729 {
730 DNDS_assert_info(false, "no valid determination" +
731 fmt::format(" [{} {}], ", exists_single, exists_rank_offsets) + name);
732 }
733 }
734 else
735 {
736 DNDS_assert_info(false, "no valid determination; arrays of ArrayGlobalOffset_One need to be explicitly designated" +
737 fmt::format(" [{} {}], ", exists_single, exists_rank_offsets) + name);
738 }
739 }
740 else
741 DNDS_assert_info(false, "no valid determination; no dataset found " +
742 fmt::format(" [{} {}], ", exists_single, exists_rank_offsets) + name);
743 }
744 else
745 {
746 if (offset == ArrayGlobalOffset_One)
747 {
748 index nGlobal{-1};
749 int dim2{-1};
750 if (v != nullptr)
751 {
752 H5_ReadDataset(group_id, name.c_str(), nGlobal, 0, size, T_H5TYPE, dxpl_id, dapl_id, v, dim2);
753 DNDS_assert(dim2 == -1);
754 DNDS_assert(size == nGlobal);
755 }
756 else
757 {
758 H5_ReadDataset(group_id, name.c_str(), nGlobal, -1, -1, T_H5TYPE, dxpl_id, dapl_id, nullptr, dim2);
759 DNDS_assert(dim2 == -1);
760 size = nGlobal;
761 }
762 }
763 // else if (offset == ArrayGlobalOffset_Parts)
764 // {
765 // index nGlobal{-1};
766 // int dim2{-1};
767 // if (v != nullptr)
768 // {
769 // H5_ReadDataset(group_id, (name + "::" + std::to_string(mpi.rank)).c_str(), nGlobal, 0, size, T_H5TYPE, dxpl_id, dapl_id, v, dim2);
770 // DNDS_assert(dim2 == -1);
771 // DNDS_assert_info(size == nGlobal, fmt::format(", {}, {}", size, nGlobal));
772 // }
773 // else
774 // {
775 // H5_ReadDataset(group_id, (name + "::" + std::to_string(mpi.rank)).c_str(), nGlobal, -1, -1, T_H5TYPE, dxpl_id, dapl_id, nullptr, dim2);
776 // DNDS_assert(dim2 == -1);
777 // size = nGlobal;
778 // }
779 // }
780 else if (offset == ArrayGlobalOffset_EvenSplit)
781 {
782 // Even-split read: each rank reads ~N_global/nRanks rows.
783 // First pass (v==nullptr): query global size, compute local range, set offset.
784 // Second pass (v!=nullptr): read the data slice.
785 index nGlobal{-1};
786 int dim2{-1};
787 if (v == nullptr)
788 {
789 H5_ReadDataset(group_id, name.c_str(), nGlobal, -1, -1, T_H5TYPE, dxpl_id, dapl_id, nullptr, dim2);
790 DNDS_assert(dim2 == -1);
791 auto [start, end] = EvenSplitRange(mpi.rank, mpi.size, nGlobal);
792 size = end - start;
793 offset = ArrayGlobalOffset{index(size), start};
794 }
795 else
796 {
797 DNDS_assert_info(offset.isDist(), "EvenSplit second pass must have a resolved offset");
798 H5_ReadDataset(group_id, name.c_str(), nGlobal, offset.offset(), size, T_H5TYPE, dxpl_id, dapl_id, v, dim2);
799 DNDS_assert(dim2 == -1);
800 }
801 }
802 else if (offset.isDist())
803 {
804 if (v == nullptr)
805 {
806 // First pass (size query): local size is already in offset.size().
807 size = offset.size();
808 }
809 else
810 {
811 // Second pass (data read): read `size` elements at offset.offset().
812 // size==0 is valid for ranks with empty partitions (nGlobal < nRanks).
813 index nGlobal{-1};
814 int dim2{-1};
815 H5_ReadDataset(group_id, name.c_str(), nGlobal, offset.offset(), size, T_H5TYPE, dxpl_id, dapl_id, v, dim2);
816 DNDS_assert(dim2 == -1);
817 }
818 }
819 else
820 DNDS_assert_info(false, "offset ill-formed");
821 }
822
823 herr = H5Pclose(dxpl_id), H5CHECK_Close;
824 herr = H5Pclose(dapl_id), H5CHECK_Close;
825 herr = H5Gclose(group_id), H5CHECK_Close;
826 }
827
828 void SerializerH5::ReadIndexVector(const std::string &name, std::vector<index> &v, ArrayGlobalOffset &offset)
829 {
830 size_t size{0};
831 index dummy{};
832 ReadDataVector<index>(name, nullptr, size, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
833 v.resize(size);
834 DNDS_assert(!(offset == ArrayGlobalOffset_Unknown));
835 ReadDataVector<index>(name, size == 0 ? &dummy : v.data(), size, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
836 }
837 void SerializerH5::ReadRowsizeVector(const std::string &name, std::vector<rowsize> &v, ArrayGlobalOffset &offset)
838 {
839 size_t size{0};
840 rowsize dummy{};
841 ReadDataVector<rowsize>(name, nullptr, size, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
842 v.resize(size);
843 DNDS_assert(!(offset == ArrayGlobalOffset_Unknown));
844 ReadDataVector<rowsize>(name, size == 0 ? &dummy : v.data(), size, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
845 }
846 void SerializerH5::ReadRealVector(const std::string &name, std::vector<real> &v, ArrayGlobalOffset &offset)
847 {
848 size_t size{0};
849 real dummy{};
850 ReadDataVector<real>(name, nullptr, size, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
851 v.resize(size);
852 DNDS_assert(!(offset == ArrayGlobalOffset_Unknown));
853 ReadDataVector<real>(name, size == 0 ? &dummy : v.data(), size, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
854 }
855 void SerializerH5::ReadString(const std::string &name, std::string &v)
856 {
857 ReadAttributeScalar<std::string>(name, v, h5file, reading, cP, collectiveMetadataRW);
858 }
860 {
861 using tValue = host_device_vector<index>;
862 herr_t herr;
863 std::string refPath;
864 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, collectiveMetadataRW);
865 htri_t exists_ref = H5Aexists(group_id, (name + "::ref").c_str());
866 herr = H5Gclose(group_id), H5CHECK_Close;
867 if (exists_ref > 0)
868 {
869 this->ReadString(name + "::ref", refPath);
870 }
871 else
872 {
873 refPath = cP + "/" + name;
874 }
875
876 if (pth_2_ssp.count(refPath))
877 {
878 v = *((ssp<tValue> *)(pth_2_ssp[refPath])); // ! reform this (and in json counterpart) to use reinterpret_cast or use STL's tools
879 }
880 else
881 {
882 v = std::make_shared<tValue>();
883 pth_2_ssp[refPath] = &v;
884
885 size_t size;
886 index dummy{};
887 ReadDataVector<index>(refPath, nullptr, size, offset, h5file, reading, "/", mpi, collectiveMetadataRW, collectiveDataRW);
888 v->resize(size);
889 DNDS_assert(!(offset == ArrayGlobalOffset_Unknown));
890 ReadDataVector<index>(refPath, size == 0 ? &dummy : v->data(), size, offset, h5file, reading, "/", mpi, collectiveMetadataRW, collectiveDataRW);
891 }
892 }
894 {
895 using tValue = host_device_vector<rowsize>;
896 herr_t herr;
897 std::string refPath;
898 hid_t group_id = GetGroupOfFileIfExist(h5file, reading, cP, collectiveMetadataRW);
899 htri_t exists_ref = H5Aexists(group_id, (name + "::ref").c_str());
900 herr = H5Gclose(group_id), H5CHECK_Close;
901 if (exists_ref > 0)
902 {
903 this->ReadString(name + "::ref", refPath);
904 }
905 else
906 {
907 refPath = cP + "/" + name;
908 }
909
910 if (pth_2_ssp.count(refPath))
911 {
912 v = *((ssp<tValue> *)(pth_2_ssp[refPath])); // ! reform this (and in json counterpart) to use reinterpret_cast or use STL's tools
913 }
914 else
915 {
916 v = std::make_shared<tValue>();
917 pth_2_ssp[refPath] = &v;
918
919 size_t size;
920 rowsize dummy{};
921 ReadDataVector<rowsize>(refPath, nullptr, size, offset, h5file, reading, "/", mpi, collectiveMetadataRW, collectiveDataRW);
922 v->resize(size);
923 DNDS_assert(!(offset == ArrayGlobalOffset_Unknown));
924 ReadDataVector<rowsize>(refPath, size == 0 ? &dummy : v->data(), size, offset, h5file, reading, "/", mpi, collectiveMetadataRW, collectiveDataRW);
925 }
926 }
927 void SerializerH5::WriteUint8Array(const std::string &name, const uint8_t *data, index size, ArrayGlobalOffset offset)
928 {
929 WriteDataVector<uint8_t>(name, data, size, offset, chunksize, deflateLevel, h5file, reading, cP, mpi, commDup, collectiveMetadataRW, collectiveDataRW);
930 }
931 void SerializerH5::ReadUint8Array(const std::string &name, uint8_t *data, index &size, ArrayGlobalOffset &offset)
932 {
933 size_t size_size_t{0};
934 if (data == nullptr)
935 {
936 ReadDataVector<uint8_t>(name, nullptr, size_size_t, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
937 size = size_size_t; // todo: check overflow?
938 DNDS_assert(!(offset == ArrayGlobalOffset_Unknown));
939 }
940 else
941 {
942 DNDS_assert(!(offset == ArrayGlobalOffset_Unknown));
943 size_size_t = size;
944 ReadDataVector<uint8_t>(name, data, size_size_t, offset, h5file, reading, cP, mpi, collectiveMetadataRW, collectiveDataRW);
945 }
946 }
947}
#define DNDS_assert_info(expr, info)
Debug-only assertion with an extra std::string info message.
Definition Errors.hpp:113
#define DNDS_assert(expr)
Debug-only assertion (compiled out when DNDS_NDEBUG is defined). Prints the expression + file/line + ...
Definition Errors.hpp:108
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:150
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:220
void AllreduceOneIndex(index &v, MPI_Op op, const MPIInfo &mpi)
Single-scalar Allreduce helper for indices (in-place, count = 1).
Definition MPI.hpp:679
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:109
std::pair< index, index > EvenSplitRange(int rank, int nRanks, index nGlobal)
Split a global range [0, nGlobal) evenly among nRanks workers.
Definition Defines.hpp:129
std::vector< std::string > splitSString(const std::string &str, char delim)
Definition Defines.hpp:772
int64_t index
Global row / DOF index type (signed 64-bit; handles multi-billion-cell meshes).
Definition Defines.hpp:107
std::shared_ptr< T > ssp
Shortened alias for std::shared_ptr used pervasively in DNDSR.
Definition Defines.hpp:138
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:105
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:49
Lightweight bundle of an MPI communicator and the calling rank's coordinates.
Definition MPI.hpp:215
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