DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
test_Serializer.cpp
Go to the documentation of this file.
1/**
2 * @file test_Serializer.cpp
3 * @brief Doctest-based unit tests for DNDS Serializer classes.
4 *
5 * Covers SerializerJSON (per-rank, no MPI needed for logic) and
6 * SerializerH5 (parallel HDF5, requires MPI). The JSON tests verify
7 * scalar, vector, uint8 (with/without codec), path operations, and
8 * shared-pointer deduplication. The H5 tests verify scalar, vector,
9 * distributed vector (non-uniform per-rank sizes), uint8 (two-pass read),
10 * path operations, and string round-trips.
11 *
12 * @note H5 parallel I/O requires all MPI ranks to open the same file.
13 * The TmpH5() helper broadcasts rank 0's PID so that the filename is
14 * identical across all ranks. FileGuard uses shared=true for H5 files
15 * so only rank 0 deletes the file, followed by an MPI_Barrier.
16 *
17 * Build: `cmake --build . -t dnds_test_serializer -j8`
18 * Run: `mpirun -np 1 ./dnds_test_serializer` (JSON + H5 single rank)
19 * `mpirun -np 4 ./dnds_test_serializer` (H5 distributed tests)
20 *
21 * @see @ref dnds_unit_tests for the full test-suite overview.
22 * @test SerializerJSON scalar/vector/uint8/path/pointer round-trip,
23 * SerializerH5 scalar/vector/distributed/uint8/path/string round-trip.
24 */
25
26#define DOCTEST_CONFIG_IMPLEMENT
27#include "doctest.h"
29#include "DNDS/SerializerH5.hpp"
30#include "DNDS/ArrayPair.hpp"
31#include "DNDS/MPI.hpp"
32#include <filesystem>
33#include <numeric>
34#include <algorithm>
35
36using namespace DNDS;
37namespace S = DNDS::Serializer;
38namespace fs = std::filesystem;
39
40int main(int argc, char **argv)
41{
42 MPI_Init(&argc, &argv);
43 doctest::Context ctx;
44 ctx.applyCommandLine(argc, argv);
45 int res = ctx.run();
46 MPI_Finalize();
47 return res;
48}
49
50// Helper: build a unique temp file name that won't collide across ranks or runs.
51static std::string TmpJSON(const std::string &tag)
52{
53 int rank = 0;
54 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
55 return fmt::format("__test_ser_{}_{}_r{}.json", tag, ::getpid(), rank);
56}
57
58static std::string TmpH5(const std::string &tag)
59{
60 // H5 parallel I/O requires ALL ranks to open the same file.
61 // Broadcast rank 0's PID so the filename is identical everywhere.
62 int pid = 0;
63 int rank = 0;
64 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
65 if (rank == 0)
66 pid = ::getpid();
67 MPI_Bcast(&pid, 1, MPI_INT, 0, MPI_COMM_WORLD);
68 return fmt::format("__test_ser_{}_{}.h5", tag, pid);
69}
70
71// RAII guard that removes a file on destruction.
72// For H5 files (shared across ranks), only rank 0 should remove.
74{
75 std::string path;
76 bool ownerOnly; // if true, only rank 0 removes
77 explicit FileGuard(std::string p, bool shared = false)
78 : path(std::move(p)), ownerOnly(shared)
79 {
80 }
82 {
83 if (ownerOnly)
84 {
85 int rank = 0;
86 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
87 if (rank == 0)
88 fs::remove(path);
89 MPI_Barrier(MPI_COMM_WORLD);
90 }
91 else
92 {
93 fs::remove(path);
94 }
95 }
96};
97
98// ===================================================================
99// SerializerJSON — scalar round-trip
100// ===================================================================
101TEST_CASE("SerializerJSON scalar round-trip")
102{
103 std::string fname = TmpJSON("scalar");
104 FileGuard guard(fname);
105
106 // --- Write ---
107 {
109 ser.OpenFile(fname, false);
110 ser.CreatePath("/test");
111 ser.GoToPath("/test");
112
113 ser.WriteInt("myInt", 42);
114 ser.WriteIndex("myIdx", 123456789LL);
115 ser.WriteReal("myReal", 3.14);
116 ser.WriteString("myStr", "hello");
117
118 ser.CloseFile();
119 }
120}
121
122// ===================================================================
123// ArrayPair redistribute — same np round-trip
124// ===================================================================
125TEST_CASE("ArrayPair redistribute — same np round-trip")
126{
127 MPIInfo mpi(MPI_COMM_WORLD);
128 std::string fname = TmpH5("h5_redist");
129 FileGuard guard(fname, true);
130
131 DNDS::index nLocal = 100 + mpi.rank * 13;
132 std::vector<DNDS::index> origIndex(nLocal);
133 {
134 DNDS::index globalOffset = 0;
135 for (int r = 0; r < mpi.rank; r++)
136 globalOffset += 100 + r * 13;
137 for (DNDS::index i = 0; i < nLocal; i++)
138 origIndex[i] = globalOffset + (nLocal - 1 - i);
139 }
140
141 {
142 using TArray = ParArray<DNDS::real, 3>;
143 using TPair = ArrayPair<TArray>;
144 TPair pair;
145 pair.InitPair("redist::pair", mpi);
146 pair.father->Resize(nLocal);
147 pair.son->Resize(0);
148 for (DNDS::index i = 0; i < nLocal; i++)
149 for (DNDS::rowsize j = 0; j < 3; j++)
150 pair.father->operator()(i, j) = DNDS::real(origIndex[i]) * 10.0 + DNDS::real(j);
151 auto ser = std::make_shared<S::SerializerH5>(mpi);
152 ser->OpenFile(fname, false);
153 pair.WriteSerialize(ser, "data", origIndex, false, false);
154 ser->CloseFile();
155 }
156 MPI_Barrier(MPI_COMM_WORLD);
157 {
158 using TArray = ParArray<DNDS::real, 3>;
159 using TPair = ArrayPair<TArray>;
160 TPair pair;
161 pair.InitPair("redist::readPair", mpi);
162 pair.father->Resize(nLocal);
163 pair.son->Resize(0);
164 auto ser = std::make_shared<S::SerializerH5>(mpi);
165 ser->OpenFile(fname, true);
166 pair.ReadSerializeRedistributed(ser, "data", origIndex);
167 ser->CloseFile();
168 for (DNDS::index i = 0; i < nLocal; i++)
169 for (DNDS::rowsize j = 0; j < 3; j++)
170 CHECK(pair.father->operator()(i, j) == doctest::Approx(DNDS::real(origIndex[i]) * 10.0 + DNDS::real(j)));
171 }
172}
173
174TEST_CASE("ArrayPair redistribute — shuffled partition same np")
175{
176 MPIInfo mpi(MPI_COMM_WORLD);
177 std::string fname = TmpH5("h5_redist_shuffle");
178 FileGuard guard(fname, true);
179
180 DNDS::index nLocal = 50 + mpi.rank * 7;
181 DNDS::index nGlobal = 0;
182 {
183 DNDS::index tmp = nLocal;
184 MPI::Allreduce(&tmp, &nGlobal, 1, DNDS_MPI_INDEX, MPI_SUM, mpi.comm);
185 }
186
187 std::vector<DNDS::index> writeOrigIndex(nLocal);
188 {
189 DNDS::index globalOffset = 0;
190 for (int r = 0; r < mpi.rank; r++)
191 globalOffset += 50 + r * 7;
192 for (DNDS::index i = 0; i < nLocal; i++)
193 writeOrigIndex[i] = globalOffset + i;
194 }
195
196 {
197 using TArray = ParArray<DNDS::real, 2>;
198 using TPair = ArrayPair<TArray>;
199 TPair pair;
200 pair.InitPair("redistShuf::pair", mpi);
201 pair.father->Resize(nLocal);
202 pair.son->Resize(0);
203 for (DNDS::index i = 0; i < nLocal; i++)
204 for (DNDS::rowsize j = 0; j < 2; j++)
205 pair.father->operator()(i, j) = DNDS::real(writeOrigIndex[i]) * 100.0 + DNDS::real(j);
206 auto ser = std::make_shared<S::SerializerH5>(mpi);
207 ser->OpenFile(fname, false);
208 pair.WriteSerialize(ser, "data", writeOrigIndex, false, false);
209 ser->CloseFile();
210 }
211 MPI_Barrier(MPI_COMM_WORLD);
212 {
213 std::vector<DNDS::index> readOrigIndex(nLocal);
214 {
215 DNDS::index globalEndOffset = nGlobal;
216 for (int r = 0; r < mpi.rank; r++)
217 globalEndOffset -= (50 + r * 7);
218 for (DNDS::index i = 0; i < nLocal; i++)
219 readOrigIndex[i] = globalEndOffset - nLocal + i;
220 }
221 using TArray = ParArray<DNDS::real, 2>;
222 using TPair = ArrayPair<TArray>;
223 TPair pair;
224 pair.InitPair("redistShuf::readPair", mpi);
225 pair.father->Resize(nLocal);
226 pair.son->Resize(0);
227 auto ser = std::make_shared<S::SerializerH5>(mpi);
228 ser->OpenFile(fname, true);
229 pair.ReadSerializeRedistributed(ser, "data", readOrigIndex);
230 ser->CloseFile();
231 for (DNDS::index i = 0; i < nLocal; i++)
232 for (DNDS::rowsize j = 0; j < 2; j++)
233 CHECK(pair.father->operator()(i, j) == doctest::Approx(DNDS::real(readOrigIndex[i]) * 100.0 + DNDS::real(j)));
234 }
235}
236
237// ===================================================================
238// Parametric redistribution tests (type x layout x row-size)
239// ===================================================================
240
241template <class T>
243{
244 return static_cast<T>(origIdx * 100 + col * 3 + 7);
245}
246
247static int g_redist_ctr = 0;
248static std::string NextTag() { return "rd_" + std::to_string(g_redist_ctr++); }
249
250static void BuildSeqOrig(std::vector<DNDS::index> &v, DNDS::index n, int rank,
251 std::function<DNDS::index(int)> nFor)
252{
253 v.resize(n);
254 DNDS::index off = 0;
255 for (int r = 0; r < rank; r++)
256 off += nFor(r);
257 for (DNDS::index i = 0; i < n; i++)
258 v[i] = off + i;
259}
260
261static void BuildRevOrig(std::vector<DNDS::index> &v, DNDS::index n,
262 DNDS::index nGlobal, int rank,
263 std::function<DNDS::index(int)> nFor)
264{
265 v.resize(n);
266 DNDS::index end = nGlobal;
267 for (int r = 0; r < rank; r++)
268 end -= nFor(r);
269 for (DNDS::index i = 0; i < n; i++)
270 v[i] = end - n + i;
271}
272
274{
275};
276struct LayoutDynamic
277{
278};
279struct LayoutCSR
280{
281};
282
283template <class T, class Layout, DNDS::rowsize RS>
285{
286 using type = T;
287 using layout = Layout;
288 static constexpr DNDS::rowsize rs = RS;
289};
290
291#define REDIST_TAG_STR(T, L, RS) TYPE_TO_STRING(RedistTag<T, L, RS>)
327#undef REDIST_TAG_STR
328
329#define REDIST_ALL_TAGS \
330 RedistTag<DNDS::real, LayoutStaticFixed, 1>, RedistTag<DNDS::real, LayoutStaticFixed, 3>, RedistTag<DNDS::real, LayoutStaticFixed, 7>, \
331 RedistTag<DNDS::real, LayoutDynamic, 1>, RedistTag<DNDS::real, LayoutDynamic, 3>, RedistTag<DNDS::real, LayoutDynamic, 7>, \
332 RedistTag<DNDS::real, LayoutCSR, 0>, \
333 RedistTag<DNDS::index, LayoutStaticFixed, 1>, RedistTag<DNDS::index, LayoutStaticFixed, 3>, RedistTag<DNDS::index, LayoutStaticFixed, 7>, \
334 RedistTag<DNDS::index, LayoutDynamic, 1>, RedistTag<DNDS::index, LayoutDynamic, 3>, RedistTag<DNDS::index, LayoutDynamic, 7>, \
335 RedistTag<DNDS::index, LayoutCSR, 0>, \
336 RedistTag<uint16_t, LayoutStaticFixed, 1>, RedistTag<uint16_t, LayoutStaticFixed, 3>, RedistTag<uint16_t, LayoutStaticFixed, 7>, \
337 RedistTag<uint16_t, LayoutDynamic, 1>, RedistTag<uint16_t, LayoutDynamic, 3>, RedistTag<uint16_t, LayoutDynamic, 7>, \
338 RedistTag<uint16_t, LayoutCSR, 0>, \
339 RedistTag<int32_t, LayoutStaticFixed, 1>, RedistTag<int32_t, LayoutStaticFixed, 3>, RedistTag<int32_t, LayoutStaticFixed, 7>, \
340 RedistTag<int32_t, LayoutDynamic, 1>, RedistTag<int32_t, LayoutDynamic, 3>, RedistTag<int32_t, LayoutDynamic, 7>, \
341 RedistTag<int32_t, LayoutCSR, 0>, \
342 RedistTag<uint8_t, LayoutStaticFixed, 1>, RedistTag<uint8_t, LayoutStaticFixed, 3>, RedistTag<uint8_t, LayoutStaticFixed, 7>, \
343 RedistTag<uint8_t, LayoutDynamic, 1>, RedistTag<uint8_t, LayoutDynamic, 3>, RedistTag<uint8_t, LayoutDynamic, 7>, \
344 RedistTag<uint8_t, LayoutCSR, 0>
345
347{
348 using T = typename Tag::type;
349 using L = typename Tag::layout;
350 constexpr DNDS::rowsize RS = Tag::rs;
351
352 MPIInfo mpi(MPI_COMM_WORLD);
353
354 for (DNDS::index baseSize : {5, 30, 100})
355 {
356 CAPTURE(baseSize);
357 std::string fname = TmpH5(NextTag());
358 FileGuard guard(fname, true);
359
360 auto nLocalFor = [baseSize](int r) -> DNDS::index
361 { return baseSize + r * 9; };
362 DNDS::index nLocal = nLocalFor(mpi.rank);
363 DNDS::index nGlobal = 0;
364 {
365 DNDS::index tmp = nLocal;
366 MPI::Allreduce(&tmp, &nGlobal, 1, DNDS_MPI_INDEX, MPI_SUM, mpi.comm);
367 }
368
369 auto csrRowSize = [](DNDS::index origIdx) -> DNDS::rowsize
370 { return DNDS::rowsize(1 + origIdx % 7); };
371
372 std::vector<DNDS::index> writeOrig;
373 BuildSeqOrig(writeOrig, nLocal, mpi.rank, nLocalFor);
374
375 // ---- Write ----
376 if constexpr (std::is_same_v<L, LayoutStaticFixed>)
377 {
378 using TArray = ParArray<T, RS>;
379 using TPair = ArrayPair<TArray>;
380 TPair pair;
381 pair.InitPair("rd::wPair", mpi);
382 pair.father->Resize(nLocal);
383 pair.son->Resize(0);
384 for (DNDS::index i = 0; i < nLocal; i++)
385 for (DNDS::rowsize j = 0; j < RS; j++)
386 pair.father->operator()(i, j) = MakeTestValue<T>(writeOrig[i], j);
387 auto ser = std::make_shared<S::SerializerH5>(mpi);
388 ser->OpenFile(fname, false);
389 pair.WriteSerialize(ser, "data", writeOrig, false, false);
390 ser->CloseFile();
391 }
392 else if constexpr (std::is_same_v<L, LayoutDynamic>)
393 {
394 using TArray = ParArray<T, DynamicSize>;
395 using TPair = ArrayPair<TArray>;
396 TPair pair;
397 pair.InitPair("rd::wPair", mpi);
398 pair.father->Resize(nLocal, RS);
399 pair.son->Resize(0, RS);
400 for (DNDS::index i = 0; i < nLocal; i++)
401 for (DNDS::rowsize j = 0; j < RS; j++)
402 pair.father->operator()(i, j) = MakeTestValue<T>(writeOrig[i], j);
403 auto ser = std::make_shared<S::SerializerH5>(mpi);
404 ser->OpenFile(fname, false);
405 pair.WriteSerialize(ser, "data", writeOrig, false, false);
406 ser->CloseFile();
407 }
408 else
409 {
410 using TArray = ParArray<T, NonUniformSize>;
411 using TPair = ArrayPair<TArray>;
412 TPair pair;
413 pair.InitPair("rd::wPair", mpi);
414 pair.father->Resize(nLocal);
415 pair.son->Resize(0);
416 for (DNDS::index i = 0; i < nLocal; i++)
417 pair.father->ResizeRow(i, csrRowSize(writeOrig[i]));
418 pair.father->Compress();
419 for (DNDS::index i = 0; i < nLocal; i++)
420 for (DNDS::rowsize j = 0; j < pair.father->RowSize(i); j++)
421 pair.father->operator()(i, j) = MakeTestValue<T>(writeOrig[i], j);
422 auto ser = std::make_shared<S::SerializerH5>(mpi);
423 ser->OpenFile(fname, false);
424 pair.WriteSerialize(ser, "data", writeOrig, false, false);
425 ser->CloseFile();
426 }
427 MPI_Barrier(MPI_COMM_WORLD);
428
429 // ---- Read with reversed partition ----
430 std::vector<DNDS::index> readOrig;
431 BuildRevOrig(readOrig, nLocal, nGlobal, mpi.rank, nLocalFor);
432
433 if constexpr (std::is_same_v<L, LayoutStaticFixed>)
434 {
435 using TArray = ParArray<T, RS>;
436 using TPair = ArrayPair<TArray>;
437 TPair pair;
438 pair.InitPair("rd::rPair", mpi);
439 pair.father->Resize(nLocal);
440 pair.son->Resize(0);
441 auto ser = std::make_shared<S::SerializerH5>(mpi);
442 ser->OpenFile(fname, true);
443 pair.ReadSerializeRedistributed(ser, "data", readOrig);
444 ser->CloseFile();
445 for (DNDS::index i = 0; i < nLocal; i++)
446 for (DNDS::rowsize j = 0; j < RS; j++)
447 CHECK(pair.father->operator()(i, j) == MakeTestValue<T>(readOrig[i], j));
448 }
449 else if constexpr (std::is_same_v<L, LayoutDynamic>)
450 {
451 using TArray = ParArray<T, DynamicSize>;
452 using TPair = ArrayPair<TArray>;
453 TPair pair;
454 pair.InitPair("rd::rPair", mpi);
455 pair.father->Resize(nLocal, RS);
456 pair.son->Resize(0, RS);
457 auto ser = std::make_shared<S::SerializerH5>(mpi);
458 ser->OpenFile(fname, true);
459 pair.ReadSerializeRedistributed(ser, "data", readOrig);
460 ser->CloseFile();
461 for (DNDS::index i = 0; i < nLocal; i++)
462 for (DNDS::rowsize j = 0; j < RS; j++)
463 CHECK(pair.father->operator()(i, j) == MakeTestValue<T>(readOrig[i], j));
464 }
465 else
466 {
467 using TArray = ParArray<T, NonUniformSize>;
468 using TPair = ArrayPair<TArray>;
469 TPair pair;
470 pair.InitPair("rd::rPair", mpi);
471 pair.father->Resize(nLocal);
472 pair.son->Resize(0);
473 auto ser = std::make_shared<S::SerializerH5>(mpi);
474 ser->OpenFile(fname, true);
475 pair.ReadSerializeRedistributed(ser, "data", readOrig);
476 ser->CloseFile();
477 for (DNDS::index i = 0; i < nLocal; i++)
478 {
479 CHECK(pair.father->RowSize(i) == csrRowSize(readOrig[i]));
480 for (DNDS::rowsize j = 0; j < csrRowSize(readOrig[i]); j++)
481 CHECK(pair.father->operator()(i, j) == MakeTestValue<T>(readOrig[i], j));
482 }
483 }
484 } // for baseSize
485}
486
487// ===================================================================
488// Different-np redistribution test (all types, layouts, sizes)
489// ===================================================================
490// Writes data using a subset of ranks (npWrite), reads back with all ranks.
491// Reuses the same REDIST_ALL_TAGS cross-product as same-np tests.
492// Requires world np >= 3. Runtime loop over npWrite and baseSize.
493
494template <class T, class L, DNDS::rowsize RS>
495void TestDifferentNp(int npWrite, DNDS::index baseSize)
496{
497 MPIInfo worldMpi(MPI_COMM_WORLD);
498 if (worldMpi.size < 3)
499 return;
500 DNDS_assert(npWrite >= 1 && npWrite < worldMpi.size);
501
502 auto csrRowSize = [](DNDS::index origIdx) -> DNDS::rowsize
503 { return DNDS::rowsize(1 + origIdx % 7); };
504
505 std::string fname = TmpH5(NextTag());
506 FileGuard guard(fname, true);
507
508 // --- Write phase: only first npWrite ranks ---
509 DNDS::index nGlobalWrite = 0;
510 {
511 int color = (worldMpi.rank < npWrite) ? 0 : MPI_UNDEFINED;
512 MPI_Comm writeComm = MPI_COMM_NULL;
513 MPI_Comm_split(MPI_COMM_WORLD, color, worldMpi.rank, &writeComm);
514
515 // Compute nGlobalWrite on all ranks
516 for (int r = 0; r < npWrite; r++)
517 nGlobalWrite += baseSize + r * 7;
518
519 if (writeComm != MPI_COMM_NULL)
520 {
521 MPIInfo writeMpi(writeComm);
522 auto nLocalFor = [baseSize](int r) -> DNDS::index
523 { return baseSize + r * 7; };
524 DNDS::index nLocal = nLocalFor(writeMpi.rank);
525
526 std::vector<DNDS::index> writeOrig;
527 BuildSeqOrig(writeOrig, nLocal, writeMpi.rank, nLocalFor);
528
529 if constexpr (std::is_same_v<L, LayoutStaticFixed>)
530 {
531 using TArray = ParArray<T, RS>;
532 using TPair = ArrayPair<TArray>;
533 TPair pair;
534 pair.InitPair("rdNp::wPair", writeMpi);
535 pair.father->Resize(nLocal);
536 pair.son->Resize(0);
537 for (DNDS::index i = 0; i < nLocal; i++)
538 for (DNDS::rowsize j = 0; j < RS; j++)
539 pair.father->operator()(i, j) = MakeTestValue<T>(writeOrig[i], j);
540 auto ser = std::make_shared<S::SerializerH5>(writeMpi);
541 ser->OpenFile(fname, false);
542 pair.WriteSerialize(ser, "data", writeOrig, false, false);
543 ser->CloseFile();
544 }
545 else if constexpr (std::is_same_v<L, LayoutDynamic>)
546 {
547 using TArray = ParArray<T, DynamicSize>;
548 using TPair = ArrayPair<TArray>;
549 TPair pair;
550 pair.InitPair("rdNp::wPair", writeMpi);
551 pair.father->Resize(nLocal, RS);
552 pair.son->Resize(0, RS);
553 for (DNDS::index i = 0; i < nLocal; i++)
554 for (DNDS::rowsize j = 0; j < RS; j++)
555 pair.father->operator()(i, j) = MakeTestValue<T>(writeOrig[i], j);
556 auto ser = std::make_shared<S::SerializerH5>(writeMpi);
557 ser->OpenFile(fname, false);
558 pair.WriteSerialize(ser, "data", writeOrig, false, false);
559 ser->CloseFile();
560 }
561 else // CSR
562 {
563 using TArray = ParArray<T, NonUniformSize>;
564 using TPair = ArrayPair<TArray>;
565 TPair pair;
566 pair.InitPair("rdNp::wPair", writeMpi);
567 pair.father->Resize(nLocal);
568 pair.son->Resize(0);
569 for (DNDS::index i = 0; i < nLocal; i++)
570 pair.father->ResizeRow(i, csrRowSize(writeOrig[i]));
571 pair.father->Compress();
572 for (DNDS::index i = 0; i < nLocal; i++)
573 for (DNDS::rowsize j = 0; j < pair.father->RowSize(i); j++)
574 pair.father->operator()(i, j) = MakeTestValue<T>(writeOrig[i], j);
575 auto ser = std::make_shared<S::SerializerH5>(writeMpi);
576 ser->OpenFile(fname, false);
577 pair.WriteSerialize(ser, "data", writeOrig, false, false);
578 ser->CloseFile();
579 }
580 MPI_Comm_free(&writeComm);
581 }
582 }
583 MPI_Barrier(MPI_COMM_WORLD);
584
585 // --- Read phase: all ranks with even-split ---
586 {
587 auto [startRow, endRow] = EvenSplitRange(worldMpi.rank, worldMpi.size, nGlobalWrite);
588 DNDS::index nLocal = endRow - startRow;
589
590 std::vector<DNDS::index> readOrig(nLocal);
591 for (DNDS::index i = 0; i < nLocal; i++)
592 readOrig[i] = startRow + i;
593
594 if constexpr (std::is_same_v<L, LayoutStaticFixed>)
595 {
596 using TArray = ParArray<T, RS>;
597 using TPair = ArrayPair<TArray>;
598 TPair pair;
599 pair.InitPair("rdNp::rPair", worldMpi);
600 pair.father->Resize(nLocal);
601 pair.son->Resize(0);
602 auto ser = std::make_shared<S::SerializerH5>(worldMpi);
603 ser->OpenFile(fname, true);
604 pair.ReadSerializeRedistributed(ser, "data", readOrig);
605 ser->CloseFile();
606 for (DNDS::index i = 0; i < nLocal; i++)
607 for (DNDS::rowsize j = 0; j < RS; j++)
608 CHECK(pair.father->operator()(i, j) == MakeTestValue<T>(readOrig[i], j));
609 }
610 else if constexpr (std::is_same_v<L, LayoutDynamic>)
611 {
612 using TArray = ParArray<T, DynamicSize>;
613 using TPair = ArrayPair<TArray>;
614 TPair pair;
615 pair.InitPair("rdNp::rPair", worldMpi);
616 pair.father->Resize(nLocal, RS);
617 pair.son->Resize(0, RS);
618 auto ser = std::make_shared<S::SerializerH5>(worldMpi);
619 ser->OpenFile(fname, true);
620 pair.ReadSerializeRedistributed(ser, "data", readOrig);
621 ser->CloseFile();
622 for (DNDS::index i = 0; i < nLocal; i++)
623 for (DNDS::rowsize j = 0; j < RS; j++)
624 CHECK(pair.father->operator()(i, j) == MakeTestValue<T>(readOrig[i], j));
625 }
626 else // CSR
627 {
628 using TArray = ParArray<T, NonUniformSize>;
629 using TPair = ArrayPair<TArray>;
630 TPair pair;
631 pair.InitPair("rdNp::rPair", worldMpi);
632 pair.father->Resize(nLocal);
633 pair.son->Resize(0);
634 auto ser = std::make_shared<S::SerializerH5>(worldMpi);
635 ser->OpenFile(fname, true);
636 pair.ReadSerializeRedistributed(ser, "data", readOrig);
637 ser->CloseFile();
638 for (DNDS::index i = 0; i < nLocal; i++)
639 {
640 CHECK(pair.father->RowSize(i) == csrRowSize(readOrig[i]));
641 for (DNDS::rowsize j = 0; j < csrRowSize(readOrig[i]); j++)
642 CHECK(pair.father->operator()(i, j) == MakeTestValue<T>(readOrig[i], j));
643 }
644 }
645 }
646}
647
648TEST_CASE_TEMPLATE("redistribute different-np", Tag, REDIST_ALL_TAGS)
649{
650 using T = typename Tag::type;
651 using L = typename Tag::layout;
652 constexpr DNDS::rowsize RS = Tag::rs;
653
654 MPIInfo mpi(MPI_COMM_WORLD);
655 if (mpi.size < 4)
656 return;
657
658 for (int npWrite : {1, 2, 3})
659 {
660 for (DNDS::index baseSize : {5, 30, 537})
661 {
662 CAPTURE(npWrite);
663 CAPTURE(baseSize);
664 TestDifferentNp<T, L, RS>(npWrite, baseSize);
665 }
666 }
667}
668
669// ===================================================================
670// SerializerJSON — vector round-trip
671// ===================================================================
672TEST_CASE("SerializerJSON vector round-trip")
673{
674 std::string fname = TmpJSON("vec");
675 FileGuard guard(fname);
676
677 const std::vector<real> wReals = {1.1, 2.2, 3.3, 4.4, 5.5};
678 const std::vector<DNDS::index> wIndices = {10, 20, 30};
679 const std::vector<rowsize> wRows = {1, 2, 3, 4};
680
681 // --- Write ---
682 {
684 ser.OpenFile(fname, false);
685 ser.CreatePath("/vectors");
686 ser.GoToPath("/vectors");
687
688 ser.WriteRealVector("reals", wReals, S::ArrayGlobalOffset_Unknown);
689 ser.WriteIndexVector("indices", wIndices, S::ArrayGlobalOffset_Unknown);
690 ser.WriteRowsizeVector("rows", wRows, S::ArrayGlobalOffset_Unknown);
691
692 ser.CloseFile();
693 }
694
695 // --- Read ---
696 {
698 ser.OpenFile(fname, true);
699 ser.GoToPath("/vectors");
700
701 std::vector<real> rReals;
702 std::vector<DNDS::index> rIndices;
703 std::vector<rowsize> rRows;
704 S::ArrayGlobalOffset off = S::ArrayGlobalOffset_Unknown;
705
706 ser.ReadRealVector("reals", rReals, off);
707 off = S::ArrayGlobalOffset_Unknown;
708 ser.ReadIndexVector("indices", rIndices, off);
709 off = S::ArrayGlobalOffset_Unknown;
710 ser.ReadRowsizeVector("rows", rRows, off);
711
712 REQUIRE(rReals.size() == wReals.size());
713 for (size_t i = 0; i < wReals.size(); ++i)
714 CHECK(rReals[i] == doctest::Approx(wReals[i]));
715
716 REQUIRE(rIndices.size() == wIndices.size());
717 for (size_t i = 0; i < wIndices.size(); ++i)
718 CHECK(rIndices[i] == wIndices[i]);
719
720 REQUIRE(rRows.size() == wRows.size());
721 for (size_t i = 0; i < wRows.size(); ++i)
722 CHECK(rRows[i] == wRows[i]);
723
724 ser.CloseFile();
725 }
726}
727
728// ===================================================================
729// SerializerJSON — uint8 array round-trip (with and without codec)
730// ===================================================================
731TEST_CASE("SerializerJSON uint8 array round-trip")
732{
733 // Prepare a known byte pattern 0..255 repeated to 512 bytes.
734 std::vector<uint8_t> pattern(512);
735 for (size_t i = 0; i < pattern.size(); ++i)
736 pattern[i] = static_cast<uint8_t>(i & 0xFF);
737
738 SUBCASE("without codec")
739 {
740 std::string fname = TmpJSON("u8_raw");
741 FileGuard guard(fname);
742
743 {
745 ser.SetUseCodecOnUint8(false);
746 ser.OpenFile(fname, false);
747 ser.CreatePath("/data");
748 ser.GoToPath("/data");
749 ser.WriteUint8Array("bytes", pattern.data(),
750 static_cast<DNDS::index>(pattern.size()),
751 S::ArrayGlobalOffset_Unknown);
752 ser.CloseFile();
753 }
754
755 {
757 ser.OpenFile(fname, true);
758 ser.GoToPath("/data");
759
760 DNDS::index sz = 0;
761 S::ArrayGlobalOffset off = S::ArrayGlobalOffset_Unknown;
762 // First call: get size only (nullptr)
763 ser.ReadUint8Array("bytes", nullptr, sz, off);
764 REQUIRE(sz == static_cast<DNDS::index>(pattern.size()));
765
766 std::vector<uint8_t> readBack(sz);
767 off = S::ArrayGlobalOffset_Unknown;
768 ser.ReadUint8Array("bytes", readBack.data(), sz, off);
769 CHECK(readBack == pattern);
770
771 ser.CloseFile();
772 }
773 }
774
775 SUBCASE("with codec (base64 + zlib)")
776 {
777 std::string fname = TmpJSON("u8_codec");
778 FileGuard guard(fname);
779
780 {
782 ser.SetUseCodecOnUint8(true);
783 ser.SetDeflateLevel(5);
784 ser.OpenFile(fname, false);
785 ser.CreatePath("/data");
786 ser.GoToPath("/data");
787 ser.WriteUint8Array("bytes", pattern.data(),
788 static_cast<DNDS::index>(pattern.size()),
789 S::ArrayGlobalOffset_Unknown);
790 ser.CloseFile();
791 }
792
793 {
795 ser.OpenFile(fname, true);
796 ser.GoToPath("/data");
797
798 DNDS::index sz = 0;
799 S::ArrayGlobalOffset off = S::ArrayGlobalOffset_Unknown;
800 ser.ReadUint8Array("bytes", nullptr, sz, off);
801 REQUIRE(sz == static_cast<DNDS::index>(pattern.size()));
802
803 std::vector<uint8_t> readBack(sz);
804 off = S::ArrayGlobalOffset_Unknown;
805 ser.ReadUint8Array("bytes", readBack.data(), sz, off);
806 CHECK(readBack == pattern);
807
808 ser.CloseFile();
809 }
810 }
811}
812
813// ===================================================================
814// SerializerJSON — path operations
815// ===================================================================
816TEST_CASE("SerializerJSON path operations")
817{
818 std::string fname = TmpJSON("paths");
819 FileGuard guard(fname);
820
822 ser.OpenFile(fname, false);
823
824 // Create nested hierarchy
825 ser.CreatePath("/a");
826 ser.GoToPath("/a");
827 ser.CreatePath("b");
828 ser.GoToPath("b");
829 ser.CreatePath("c");
830 ser.CreatePath("d");
831
832 // Current path should be /a/b
833 CHECK(ser.GetCurrentPath() == "/a/b");
834
835 // Listing /a/b should contain "c" and "d"
836 auto entries = ser.ListCurrentPath();
837 CHECK(entries.count("c") == 1);
838 CHECK(entries.count("d") == 1);
839
840 // Go to /a/b/c and write something to verify it's valid
841 ser.GoToPath("c");
842 CHECK(ser.GetCurrentPath() == "/a/b/c");
843 ser.WriteInt("val", 99);
844
845 ser.CloseFile();
846
847 // Verify the written value survived
848 {
849 S::SerializerJSON reader;
850 reader.OpenFile(fname, true);
851 reader.GoToPath("/a/b/c");
852 int v = 0;
853 reader.ReadInt("val", v);
854 CHECK(v == 99);
855 reader.CloseFile();
856 }
857}
858
859// ===================================================================
860// SerializerJSON — shared pointer deduplication
861// ===================================================================
862TEST_CASE("SerializerJSON shared pointer deduplication")
863{
864 std::string fname = TmpJSON("shared");
865 FileGuard guard(fname);
866
867 // Build a shared index vector
868 auto sharedVec = std::make_shared<host_device_vector<DNDS::index>>();
869 sharedVec->resize(5);
870 for (size_t i = 0; i < 5; ++i)
871 (*sharedVec)[i] = static_cast<DNDS::index>(100 + i);
872
873 // --- Write: same shared_ptr under two names ---
874 {
876 ser.OpenFile(fname, false);
877 ser.CreatePath("/dedup");
878 ser.GoToPath("/dedup");
879
880 ser.WriteSharedIndexVector("first", sharedVec, S::ArrayGlobalOffset_Unknown);
881 ser.WriteSharedIndexVector("second", sharedVec, S::ArrayGlobalOffset_Unknown);
882
883 ser.CloseFile();
884 }
885
886 // --- Read back both ---
887 {
889 ser.OpenFile(fname, true);
890 ser.GoToPath("/dedup");
891
894 S::ArrayGlobalOffset off = S::ArrayGlobalOffset_Unknown;
895
896 ser.ReadSharedIndexVector("first", readFirst, off);
897 off = S::ArrayGlobalOffset_Unknown;
898 ser.ReadSharedIndexVector("second", readSecond, off);
899
900 REQUIRE(readFirst);
901 REQUIRE(readSecond);
902
903 // Both should resolve to the same underlying data (deduplication)
904 CHECK(readFirst.get() == readSecond.get());
905
906 // Values should match
907 REQUIRE(readFirst->size() == 5);
908 for (size_t i = 0; i < 5; ++i)
909 CHECK((*readFirst)[i] == static_cast<DNDS::index>(100 + i));
910
911 ser.CloseFile();
912 }
913}
914
915// ===================================================================
916// SerializerH5 — scalar round-trip
917// ===================================================================
918TEST_CASE("SerializerH5 scalar round-trip")
919{
920 MPIInfo mpi(MPI_COMM_WORLD);
921 std::string fname = TmpH5("h5scalar");
922 FileGuard guard(fname, true);
923
924 // --- Write ---
925 {
926 S::SerializerH5 ser(mpi);
927 ser.OpenFile(fname, false);
928 ser.CreatePath("/scalars");
929 ser.GoToPath("/scalars");
930
931 ser.WriteInt("myInt", 42);
932 ser.WriteIndex("myIdx", 987654321LL);
933 ser.WriteReal("myReal", 2.718281828);
934 ser.WriteString("myStr", "world");
935
936 ser.CloseFile();
937 }
938
939 MPI_Barrier(MPI_COMM_WORLD);
940
941 // --- Read ---
942 {
943 S::SerializerH5 ser(mpi);
944 ser.OpenFile(fname, true);
945 ser.GoToPath("/scalars");
946
947 int rInt = 0;
948 DNDS::index rIdx = 0;
949 real rReal = 0.0;
950 std::string rStr;
951
952 ser.ReadInt("myInt", rInt);
953 ser.ReadIndex("myIdx", rIdx);
954 ser.ReadReal("myReal", rReal);
955 ser.ReadString("myStr", rStr);
956
957 CHECK(rInt == 42);
958 CHECK(rIdx == 987654321LL);
959 CHECK(rReal == doctest::Approx(2.718281828));
960 CHECK(rStr == "world");
961
962 ser.CloseFile();
963 }
964}
965
966// ===================================================================
967// SerializerH5 — vector round-trip (each rank writes its own portion)
968// ===================================================================
969TEST_CASE("SerializerH5 vector round-trip")
970{
971 MPIInfo mpi(MPI_COMM_WORLD);
972 std::string fname = TmpH5("h5vec");
973 FileGuard guard(fname, true);
974
975 const DNDS::index N = 10; // elements per rank
976 DNDS::index globalSize = static_cast<DNDS::index>(mpi.size) * N;
977 DNDS::index myOffset = static_cast<DNDS::index>(mpi.rank) * N;
978
979 // Build local data
980 std::vector<real> wReals(N);
981 std::vector<DNDS::index> wIndices(N);
982 for (DNDS::index i = 0; i < N; ++i)
983 {
984 wReals[i] = static_cast<real>(myOffset + i) * 0.1;
985 wIndices[i] = myOffset + i;
986 }
987
988 S::ArrayGlobalOffset distOff(globalSize, myOffset);
989
990 // --- Write ---
991 {
992 S::SerializerH5 ser(mpi);
993 ser.OpenFile(fname, false);
994 ser.CreatePath("/vecs");
995 ser.GoToPath("/vecs");
996
997 ser.WriteRealVector("reals", wReals, distOff);
998 ser.WriteIndexVector("indices", wIndices, distOff);
999
1000 ser.CloseFile();
1001 }
1002
1003 MPI_Barrier(MPI_COMM_WORLD);
1004
1005 // --- Read (always use ArrayGlobalOffset_Unknown; H5 auto-detects from ::rank_offsets) ---
1006 {
1007 S::SerializerH5 ser(mpi);
1008 ser.OpenFile(fname, true);
1009 ser.GoToPath("/vecs");
1010
1011 std::vector<real> rReals;
1012 std::vector<DNDS::index> rIndices;
1013 S::ArrayGlobalOffset offR = S::ArrayGlobalOffset_Unknown;
1014 S::ArrayGlobalOffset offI = S::ArrayGlobalOffset_Unknown;
1015
1016 ser.ReadRealVector("reals", rReals, offR);
1017 ser.ReadIndexVector("indices", rIndices, offI);
1018
1019 REQUIRE(rReals.size() == static_cast<size_t>(N));
1020 for (DNDS::index i = 0; i < N; ++i)
1021 CHECK(rReals[i] == doctest::Approx(wReals[i]));
1022
1023 REQUIRE(rIndices.size() == static_cast<size_t>(N));
1024 for (DNDS::index i = 0; i < N; ++i)
1025 CHECK(rIndices[i] == wIndices[i]);
1026
1027 ser.CloseFile();
1028 }
1029}
1030
1031// ===================================================================
1032// SerializerH5 — distributed vector (unknown offset on read)
1033// ===================================================================
1034TEST_CASE("SerializerH5 distributed vector")
1035{
1036 MPIInfo mpi(MPI_COMM_WORLD);
1037 std::string fname = TmpH5("h5dist");
1038 FileGuard guard(fname, true);
1039
1040 // Each rank writes a variable-size portion.
1041 // Rank r writes (r + 1) * 3 elements.
1042 const DNDS::index localN = static_cast<DNDS::index>(mpi.rank + 1) * 3;
1043
1044 // Compute global size and per-rank offset via MPI scan.
1045 DNDS::index localN64 = localN;
1046 DNDS::index myOffset = 0;
1047 MPI_Scan(&localN64, &myOffset, 1, MPI_INT64_T, MPI_SUM, MPI_COMM_WORLD);
1048 myOffset -= localN64; // exclusive prefix sum
1050 MPI_Allreduce(&localN64, &globalSize, 1, MPI_INT64_T, MPI_SUM, MPI_COMM_WORLD);
1051
1052 S::ArrayGlobalOffset distOff(globalSize, myOffset);
1053
1054 // Fill local data: each element is (globalIndex * 7 + 3)
1055 std::vector<DNDS::index> wData(localN);
1056 for (DNDS::index i = 0; i < localN; ++i)
1057 wData[i] = (myOffset + i) * 7 + 3;
1058
1059 // --- Write ---
1060 {
1061 S::SerializerH5 ser(mpi);
1062 ser.OpenFile(fname, false);
1063 ser.CreatePath("/dist");
1064 ser.GoToPath("/dist");
1065
1066 ser.WriteIndexVector("data", wData, distOff);
1067
1068 ser.CloseFile();
1069 }
1070
1071 MPI_Barrier(MPI_COMM_WORLD);
1072
1073 // --- Read with Unknown offset (auto-detect from rank_offsets) ---
1074 {
1075 S::SerializerH5 ser(mpi);
1076 ser.OpenFile(fname, true);
1077 ser.GoToPath("/dist");
1078
1079 std::vector<DNDS::index> rData;
1080 S::ArrayGlobalOffset offR = S::ArrayGlobalOffset_Unknown;
1081
1082 ser.ReadIndexVector("data", rData, offR);
1083
1084 REQUIRE(rData.size() == static_cast<size_t>(localN));
1085 for (DNDS::index i = 0; i < localN; ++i)
1086 CHECK(rData[i] == (myOffset + i) * 7 + 3);
1087
1088 ser.CloseFile();
1089 }
1090
1091 // --- Read with explicit known offset (pre-allocate the vector) ---
1092 {
1093 S::SerializerH5 ser(mpi);
1094 ser.OpenFile(fname, true);
1095 ser.GoToPath("/dist");
1096
1097 std::vector<DNDS::index> rData;
1098 S::ArrayGlobalOffset offR = S::ArrayGlobalOffset_Unknown;
1099
1100 ser.ReadIndexVector("data", rData, offR);
1101
1102 REQUIRE(rData.size() == static_cast<size_t>(localN));
1103 for (DNDS::index i = 0; i < localN; ++i)
1104 CHECK(rData[i] == (myOffset + i) * 7 + 3);
1105
1106 ser.CloseFile();
1107 }
1108}
1109
1110// ===================================================================
1111// SerializerH5 — uint8 distributed round-trip
1112// ===================================================================
1113TEST_CASE("SerializerH5 uint8 distributed round-trip")
1114{
1115 MPIInfo mpi(MPI_COMM_WORLD);
1116 std::string fname = TmpH5("h5u8");
1117 FileGuard guard(fname, true);
1118
1119 const DNDS::index localN = 64;
1120 DNDS::index globalSize = static_cast<DNDS::index>(mpi.size) * localN;
1121 DNDS::index myOffset = static_cast<DNDS::index>(mpi.rank) * localN;
1122
1123 S::ArrayGlobalOffset distOff(globalSize, myOffset);
1124
1125 // Fill local bytes: cyclic pattern seeded by rank
1126 std::vector<uint8_t> wBytes(localN);
1127 for (DNDS::index i = 0; i < localN; ++i)
1128 wBytes[i] = static_cast<uint8_t>((myOffset + i) & 0xFF);
1129
1130 // --- Write ---
1131 {
1132 S::SerializerH5 ser(mpi);
1133 ser.OpenFile(fname, false);
1134 ser.CreatePath("/bytes");
1135 ser.GoToPath("/bytes");
1136
1137 ser.WriteUint8Array("raw", wBytes.data(), localN, distOff);
1138
1139 ser.CloseFile();
1140 }
1141
1142 MPI_Barrier(MPI_COMM_WORLD);
1143
1144 // --- Read (two-pass: nullptr to query size, then actual read with same offset) ---
1145 {
1146 S::SerializerH5 ser(mpi);
1147 ser.OpenFile(fname, true);
1148 ser.GoToPath("/bytes");
1149
1150 // Pass 1: query size (data=nullptr)
1151 DNDS::index sz = 0;
1152 S::ArrayGlobalOffset offR = S::ArrayGlobalOffset_Unknown;
1153 ser.ReadUint8Array("raw", nullptr, sz, offR);
1154 REQUIRE(sz == localN);
1155
1156 // Pass 2: actual read (reuse offR from pass 1, do NOT reset to Unknown)
1157 std::vector<uint8_t> rBytes(sz);
1158 ser.ReadUint8Array("raw", rBytes.data(), sz, offR);
1159
1160 REQUIRE(rBytes.size() == static_cast<size_t>(localN));
1161 CHECK(rBytes == wBytes);
1162
1163 ser.CloseFile();
1164 }
1165}
1166
1167// ===================================================================
1168// SerializerH5 — path and listing
1169// ===================================================================
1170TEST_CASE("SerializerH5 path operations")
1171{
1172 MPIInfo mpi(MPI_COMM_WORLD);
1173 std::string fname = TmpH5("h5paths");
1174 FileGuard guard(fname, true);
1175
1176 {
1177 S::SerializerH5 ser(mpi);
1178 ser.OpenFile(fname, false);
1179
1180 ser.CreatePath("/x");
1181 ser.GoToPath("/x");
1182
1183 CHECK(ser.GetCurrentPath() == "/x");
1184
1185 ser.WriteInt("val", 7);
1186
1187 // Groups in HDF5 are lazily created when content is written into them.
1188 // Navigate into each child and write a marker to materialize the group.
1189 auto cwd = ser.GetCurrentPath();
1190 ser.GoToPath("y");
1191 ser.WriteInt("marker", 1);
1192 ser.GoToPath(cwd);
1193 ser.GoToPath("z");
1194 ser.WriteInt("marker", 2);
1195 ser.GoToPath(cwd);
1196
1197 ser.CloseFile();
1198 }
1199
1200 MPI_Barrier(MPI_COMM_WORLD);
1201
1202 {
1203 S::SerializerH5 ser(mpi);
1204 ser.OpenFile(fname, true);
1205 ser.GoToPath("/x");
1206
1207 auto entries = ser.ListCurrentPath();
1208 CHECK(entries.count("y") == 1);
1209 CHECK(entries.count("z") == 1);
1210
1211 int v = 0;
1212 ser.ReadInt("val", v);
1213 CHECK(v == 7);
1214
1215 ser.CloseFile();
1216 }
1217}
1218
1219// ===================================================================
1220// SerializerH5 — string round-trip
1221// ===================================================================
1222TEST_CASE("SerializerH5 string round-trip")
1223{
1224 MPIInfo mpi(MPI_COMM_WORLD);
1225 std::string fname = TmpH5("h5str");
1226 FileGuard guard(fname, true);
1227
1228 {
1229 S::SerializerH5 ser(mpi);
1230 ser.OpenFile(fname, false);
1231 ser.CreatePath("/meta");
1232 ser.GoToPath("/meta");
1233 ser.WriteString("solver", "euler3D");
1234 ser.WriteString("version", "1.2.3");
1235 ser.CloseFile();
1236 }
1237
1238 MPI_Barrier(MPI_COMM_WORLD);
1239
1240 {
1241 S::SerializerH5 ser(mpi);
1242 ser.OpenFile(fname, true);
1243 ser.GoToPath("/meta");
1244
1245 std::string solver, version;
1246 ser.ReadString("solver", solver);
1247 ser.ReadString("version", version);
1248
1249 CHECK(solver == "euler3D");
1250 CHECK(version == "1.2.3");
1251
1252 ser.CloseFile();
1253 }
1254}
Father-son array pairs with device views and ghost communication.
#define DNDS_assert(expr)
Debug-only assertion (compiled out when DNDS_NDEBUG is defined). Prints the expression + file/line + ...
Definition Errors.hpp:108
MPI wrappers: MPIInfo, collective operations, type mapping, CommStrategy.
MPI-parallel HDF5 serializer implementing the SerializerBase interface.
Per-rank JSON serializer implementing the SerializerBase interface.
MPI-aware Array: adds a communicator, rank, and global index mapping.
Describes one rank's window into a globally-distributed dataset.
MPI-parallel HDF5 serializer; all ranks collectively read/write a single .h5 file.
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 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 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 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.
Per-rank JSON file serializer; each MPI rank writes its own .json file.
void ReadIndexVector(const std::string &name, std::vector< index > &v, ArrayGlobalOffset &offset) override
void ReadSharedIndexVector(const std::string &name, ssp< host_device_vector< index > > &v, ArrayGlobalOffset &offset) override
void ReadRealVector(const std::string &name, std::vector< real > &v, ArrayGlobalOffset &offset) override
std::set< std::string > ListCurrentPath() override
Names of direct children of the current path.
std::string GetCurrentPath() override
String form of the current path.
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 OpenFile(const std::string &fName, bool read) override
Open a backing file (H5 file or JSON file depending on subclass).
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 ReadInt(const std::string &name, int &v) override
Read a scalar int into v.
void ReadUint8Array(const std::string &name, uint8_t *data, index &size, ArrayGlobalOffset &offset) override
Two-pass byte array read.
void CreatePath(const std::string &p) override
Create a sub-path (H5 group / JSON object) at the current location.
void WriteRowsizeVector(const std::string &name, const std::vector< rowsize > &v, ArrayGlobalOffset offset) override
Write a rowsize vector (collective for H5).
void WriteInt(const std::string &name, int v) override
Write a scalar int under name at the current path.
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 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...
void WriteRealVector(const std::string &name, const std::vector< real > &v, ArrayGlobalOffset offset) override
Write a real vector (collective for H5).
void WriteReal(const std::string &name, real v) override
Write a scalar real under name.
void ReadRowsizeVector(const std::string &name, std::vector< rowsize > &v, ArrayGlobalOffset &offset) override
void WriteIndex(const std::string &name, index v) override
Write a scalar index under name.
void GoToPath(const std::string &p) override
Navigate to an existing path. Supports / -separated segments.
int main()
Definition dummy.cpp:3
the host side operators are provided as implemented
const MPI_Datatype DNDS_MPI_INDEX
MPI datatype matching index (= MPI_INT64_T).
Definition MPI.hpp:90
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
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
double real
Canonical floating-point scalar used throughout DNDSR (double precision).
Definition Defines.hpp:105
t_arrayDeviceView father
Convenience bundle of a father, son, and attached ArrayTransformer.
void InitPair(const std::string &name, Args &&...args)
Allocate both father and son arrays, forwarding all args to TArray constructor.
Lightweight bundle of an MPI communicator and the calling rank's coordinates.
Definition MPI.hpp:215
int size
Number of ranks in comm (-1 until initialised).
Definition MPI.hpp:221
int rank
This rank's 0-based index within comm (-1 until initialised).
Definition MPI.hpp:219
std::string path
FileGuard(std::string p, bool shared=false)
static constexpr DNDS::rowsize rs
Eigen::Matrix< real, 5, 1 > v
constexpr DNDS::index N
constexpr DNDS::index nLocal
for(DNDS::index i=0;i< nLocal;i++) for(DNDS j< dynCols;j++)(*father)(i, j)=(gOff+i) *0.1+j *0.001;DNDS::index gSize=father-> globalSize()
tVec r(NCells)
CHECK(result.size()==3)
auto res
Definition test_ODE.cpp:196
TEST_CASE("3D: VFV P2 HQM error < P1 on sinCos3D")
Eigen::Vector3d n(1.0, 0.0, 0.0)
#define REDIST_TAG_STR(T, L, RS)
#define REDIST_ALL_TAGS
TEST_CASE_TEMPLATE("redistribute", Tag, REDIST_ALL_TAGS)
void TestDifferentNp(int npWrite, DNDS::index baseSize)
T MakeTestValue(DNDS::index origIdx, DNDS::rowsize col)