DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
SerializerJSON.hpp
Go to the documentation of this file.
1#pragma once
2/// @file SerializerJSON.hpp
3/// @brief Per-rank JSON serializer implementing the SerializerBase interface.
4/// @par Unit Test Coverage (test_Serializer.cpp, MPI np=1,2,4)
5/// - Scalar round-trip: WriteInt/ReadInt, WriteIndex/ReadIndex,
6/// WriteReal/ReadReal, WriteString/ReadString
7/// - Vector round-trip: WriteRealVector, WriteIndexVector, WriteRowsizeVector
8/// - uint8 array: with and without codec (base64 + zlib)
9/// - Path operations: CreatePath, GoToPath, GetCurrentPath, ListCurrentPath
10/// - Shared pointer deduplication: WriteSharedIndexVector / ReadSharedIndexVector
11/// @par Not Yet Tested
12/// - WriteSharedRowsizeVector / ReadSharedRowsizeVector
13/// - WriteIndexVectorPerRank
14/// - SetDeflateLevel impact verification
15/// - Error handling (nonexistent file, duplicate paths)
16#include "SerializerBase.hpp"
17
18#include "DNDS/JsonUtil.hpp"
19#include <fstream>
20#include <map>
21
22namespace DNDS::Serializer
23{
24 /// @brief Per-rank JSON file serializer; each MPI rank writes its own .json file.
26 {
27 std::fstream fileStream;
28 nlohmann::json jObj;
29 bool reading = true;
30 std::vector<std::string> cPathSplit;
31 std::string cP; // current path
32 // ptr_2_pth and pth_2_ssp are inherited from SerializerBase
33
34 bool useCodecOnUint8{false};
35 int deflateLevel{5};
36
37 MPIInfo mpi; // NULL
38
39 public:
40 void SetUseCodecOnUint8(bool v) { useCodecOnUint8 = v; }
41 void SetDeflateLevel(int v) { deflateLevel = v; }
42
43 void OpenFile(const std::string &fName, bool read) override;
44 void CloseFile() override;
46 void CreatePath(const std::string &p) override;
47 void GoToPath(const std::string &p) override;
48 bool IsPerRank() override { return true; }
49 std::string GetCurrentPath() override;
50 std::set<std::string> ListCurrentPath() override;
51 int GetMPIRank() override { return -1; }
52 int GetMPISize() override { return -1; }
53 const MPIInfo &getMPI() override { return mpi; }
54
55 void WriteInt(const std::string &name, int v) override;
56 void WriteIndex(const std::string &name, index v) override;
57 void WriteReal(const std::string &name, real v) override;
58 void WriteString(const std::string &name, const std::string &v) override;
59
60 void WriteIndexVector(const std::string &name, const std::vector<index> &v, ArrayGlobalOffset offset) override;
61 void WriteRowsizeVector(const std::string &name, const std::vector<rowsize> &v, ArrayGlobalOffset offset) override;
62 void WriteRealVector(const std::string &name, const std::vector<real> &v, ArrayGlobalOffset offset) override;
63 void WriteSharedIndexVector(const std::string &name, const ssp<host_device_vector<index>> &v, ArrayGlobalOffset offset) override;
64 void WriteSharedRowsizeVector(const std::string &name, const ssp<host_device_vector<rowsize>> &v, ArrayGlobalOffset offset) override;
65
66 void WriteUint8Array(const std::string &name, const uint8_t *data, index size, ArrayGlobalOffset offset) override;
67
68 void WriteIndexVectorPerRank(const std::string &name, const std::vector<index> &v) override
69 {
70 this->WriteIndexVector(name, v, ArrayGlobalOffset_Unknown);
71 }
72
73 void ReadInt(const std::string &name, int &v) override;
74 void ReadIndex(const std::string &name, index &v) override;
75 void ReadReal(const std::string &name, real &v) override;
76 void ReadString(const std::string &name, std::string &v) override;
77
78 void ReadIndexVector(const std::string &name, std::vector<index> &v, ArrayGlobalOffset &offset) override;
79 void ReadRowsizeVector(const std::string &name, std::vector<rowsize> &v, ArrayGlobalOffset &offset) override;
80 void ReadRealVector(const std::string &name, std::vector<real> &v, ArrayGlobalOffset &offset) override;
81 void ReadSharedIndexVector(const std::string &name, ssp<host_device_vector<index>> &v, ArrayGlobalOffset &offset) override;
82 void ReadSharedRowsizeVector(const std::string &name, ssp<host_device_vector<rowsize>> &v, ArrayGlobalOffset &offset) override;
83
84 void ReadUint8Array(const std::string &name, uint8_t *data, index &size, ArrayGlobalOffset &offset) override;
85
86 ~SerializerJSON() override
87 {
89 }
90 };
91}
JSON-to-Eigen conversion utilities and nlohmann_json helper macros.
Base types and abstract interface for array serialization.
Describes one rank's window into a globally-distributed dataset.
Abstract interface for reading/writing scalars, vectors, and byte arrays.
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 ReadString(const std::string &name, std::string &v) override
Read a UTF-8 string into v.
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 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 WriteString(const std::string &name, const std::string &v) override
Write a UTF-8 string under name.
void WriteIndexVectorPerRank(const std::string &name, const std::vector< index > &v) override
Write a per-rank index vector (replicated name, independent values).
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...
int GetMPIRank() override
Rank index cached by the serializer (relevant for collective I/O).
void ReadSharedRowsizeVector(const std::string &name, ssp< host_device_vector< rowsize > > &v, ArrayGlobalOffset &offset) override
bool IsPerRank() override
Whether this serializer is per-rank (JSON-file-per-rank) or collective (shared H5 file)....
void ReadIndex(const std::string &name, index &v) override
Read a scalar index into v.
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.
const MPIInfo & getMPI() override
MPI context the serializer was opened with.
int GetMPISize() override
Rank count cached by the serializer.
void GoToPath(const std::string &p) override
Navigate to an existing path. Supports / -separated segments.
void ReadReal(const std::string &name, real &v) override
Read a scalar real into v.
Host + optional device vector of trivially copyable T.
Definition Vector.hpp:217
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
Lightweight bundle of an MPI communicator and the calling rank's coordinates.
Definition MPI.hpp:215
Eigen::Matrix< real, 5, 1 > v