MAIA bb96820c
Multiphysics at AIA
Loading...
Searching...
No Matches
parallelio_pnetcdf.h
Go to the documentation of this file.
1// Copyright (C) 2024 The m-AIA AUTHORS
2//
3// This file is part of m-AIA (https://git.rwth-aachen.de/aia/m-AIA/m-AIA)
4//
5// SPDX-License-Identifier: LGPL-3.0-only
6
7#ifndef PARALLELIO_PNETCDF_H
8#define PARALLELIO_PNETCDF_H
9
10#include <map>
11
12class ParallelIoPNetcdf : public ParallelIoBase<ParallelIoPNetcdf> {
14 // Public member functions
15 public:
16 // Constructor
17 ParallelIoPNetcdf(const MString& fileName, MInt fileMode, const MPI_Comm& mpiComm);
18 // Destructor
19 ~ParallelIoPNetcdf() override;
20 // Close the file
21 void close();
22
23 // Parallel Netcdf-specific methods (private)
24 // -> only needed inside the class (rationale: the actual data I/O calls
25 // should be hidden from the user)
26 private:
27 // Static file system-related methods
28 static MBool b_isValidFile(const MString& name, const MPI_Comm& mpiComm);
29 static const MString b_fileExt();
30
31 // File methods
32 void b_ncEndDef();
33 void b_ncRedef();
34 void b_saveHeader();
37
38 // Define mode methods
39 void b_defineArray(maiabd_type type, const MString& name, size_type noDims, size_type* totalCount);
40 void b_defineArray(maiabd_type, const MString&, const MString&, const size_type, const size_type*);
41 void b_defineScalar(maiabd_type type, const MString& name);
42
43 // Inquiry methods
44 MBool b_hasDataset(const MString& name, const size_type dimension);
45 MBool b_hasDataset(const MString& name, const MString& path);
46 MBool b_hasObject(const MString& path);
47 MInt b_getDatasetType(const MString& name);
48 void b_getDatasetNames(std::vector<MString>& names, const size_type dimension);
49 void b_getDatasetNames(std::vector<MString>& names, const MString& path);
50 void b_getGroupNames(std::vector<MString>& names, const MString& path);
52 size_type b_getDatasetNoDims(const MString& path, const MString& name);
53 size_type b_getDatasetSize(const MString& name, const size_type dimensionId = 0);
54 void b_getDatasetSize(const MString& name, const MString& path, size_type noDims, size_type* data);
55 MBool b_hasAttribute(const MString& name, const MString& datasetName = "");
56 MInt b_getAttributeType(const MString& name, const MString& datasetName = "");
57 void b_printFileHints();
58
59 // Data mode methods
60 template <class T>
61 void b_writeArray(const T* const array, const MString& name, const size_type noDims, MPI_Offset* start,
62 MPI_Offset* count, MPI_Offset memoryStride, const size_type noChunks, MPI_Offset diskStride);
63 template <class T>
64 void b_writeArray(const T*, const MString&, const MString&, const size_type, const size_type*, const size_type*,
65 const size_type*);
66 template <class T>
67 void b_writeScalar(const T scalar, const MString& name);
68 template <class T>
69 void b_readArray(T* array, const MString& name, const size_type noDims, MPI_Offset* start, MPI_Offset* count,
70 MPI_Offset memoryStride, const size_type noChunks, MPI_Offset diskStride);
71 template <class T>
72 void b_readArray(T*, const MString&, const MString&, const size_type, const size_type*, const size_type*);
73 template <class T>
74 void b_readScalar(T* scalar, const MString& name);
75
76 // Attribute methods
77 template <class T>
78 void b_setAttribute(const T* value, const MString& name, const MString& datasetName = "",
79 const size_type totalCount = 1);
80 template <class T>
81 void b_getAttribute(T* const value, const MString& name, const MString& datasetName = "",
82 const size_type totalCount = 1);
83 void b_getAttributeCount(const MString& name, size_type* totalCount, const MString& datasetName = "");
84
85 // Auxiliary methods
86 // b_error is static because it is used in b_isValidFile (also static)
87 static void b_error(MInt status, const MString& name, const MString& location);
88
89 // Parallel Netcdf-specific members
90 private:
93 struct NcDimProxy {
95 };
96 using NcDimMap = std::map<size_type, NcDimProxy>;
98};
99#endif /* PARALLELIO_H */
This class is intended to do all the heavy lifting when it comes to reading and writing "big data" fi...
Definition: parallelio.h:101
MLong size_type
Type used for all size- and offset-related values.
Definition: parallelio.h:123
maia::parallel_io::maiabd_type maiabd_type
Type used for all constants that are defined in maia::parallel_io.
Definition: parallelio.h:128
MInt b_getDatasetType(const MString &name)
Returns the data type of a dataset in the file (can be array, multi-D array or scalar).
static const MString b_fileExt()
Returns backend-specific ending of filename (either ".Netcdf" or ".Hdf5")
void b_getDatasetNames(std::vector< MString > &names, const size_type dimension)
Returns a vector of all existing datasets with the given number of dimensions in the file (if any)....
void b_writeArray(const T *const array, const MString &name, const size_type noDims, MPI_Offset *start, MPI_Offset *count, MPI_Offset memoryStride, const size_type noChunks, MPI_Offset diskStride)
Writes array data to file (generic version). [MPI]
void b_saveHeader()
Adds all additional header information that are needed in the file.
void b_defineArray(maiabd_type, const MString &, const MString &, const size_type, const size_type *)
~ParallelIoPNetcdf() override
Calls close(). [MPI]
void b_readScalar(T *scalar, const MString &name)
Read scalar data from file (generic version). [MPI]
void b_printFileHints()
Print PNetCDF file hints to cerr.
void b_writeAdditionalData()
Write additional data to file (NetCDF-specific). [MPI]
void b_getAttribute(T *const value, const MString &name, const MString &datasetName="", const size_type totalCount=1)
Retrieve an attribute from file (generic version).
std::map< size_type, NcDimProxy > NcDimMap
void b_getAttributeCount(const MString &name, size_type *totalCount, const MString &datasetName="")
void b_getGroupNames(std::vector< MString > &names, const MString &path)
static MBool b_isValidFile(const MString &name, const MPI_Comm &mpiComm)
Check if specified file is a valid HDF5 file (i.e. can be opened). [MPI]
void b_defineScalar(maiabd_type type, const MString &name)
Defines a scalar in the file. [MPI]
void b_defineArray(maiabd_type type, const MString &name, size_type noDims, size_type *totalCount)
Defines an array in the file. [MPI]
void b_setAttribute(const T *value, const MString &name, const MString &datasetName="", const size_type totalCount=1)
Set an attribute in the file (generic version).
size_type b_getDatasetNoDims(const MString &name)
Get number of dimensions of a dataset with the given name.
size_type b_getDatasetSize(const MString &name, const size_type dimensionId=0)
Get the length of one dimension of an arbitrary array in the file.
MBool b_hasObject(const MString &path)
MBool b_hasDataset(const MString &name, const size_type dimension)
Check if the file contains a dataset with the given name and number of dimensions.
MInt b_getAttributeType(const MString &name, const MString &datasetName="")
Returns the data type of an attribute in the file.
void b_addAdditionalHeader()
Write additional headers to file (e.g. grid file name, creation date etc.). [MPI]
void b_readArray(T *array, const MString &name, const size_type noDims, MPI_Offset *start, MPI_Offset *count, MPI_Offset memoryStride, const size_type noChunks, MPI_Offset diskStride)
Read array data from file (generic version). [MPI]
static void b_error(MInt status, const MString &name, const MString &location)
Check the status code of a HDF5 operation and output a meaningful message.
void close()
Close the file (normally called by the destructor but needs to be explicitely called earlier in speci...
MBool b_hasAttribute(const MString &name, const MString &datasetName="")
Check if a given attribute exists in the file.
void b_ncEndDef()
Leave the define mode (NetCDF-specific). [MPI]
void b_writeScalar(const T scalar, const MString &name)
Writes scalar data to file (generic version). [MPI]
void b_ncRedef()
Enter define mode (NetCDF-specific). [MPI]
const MString & location
Definition: functions.h:37
int32_t MInt
Definition: maiatypes.h:62
std::basic_string< char > MString
Definition: maiatypes.h:55
bool MBool
Definition: maiatypes.h:58