MAIA bb96820c
Multiphysics at AIA
Loading...
Searching...
No Matches
parallelio_hdf5.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_HDF5_H
8#define PARALLELIO_HDF5_H
9
10#ifdef MAIA_GCC_COMPILER
11#pragma GCC diagnostic ignored "-Wredundant-tags"
12#endif
13#include "hdf5.h"
14#ifdef MAIA_GCC_COMPILER
15#pragma GCC diagnostic pop
16#endif
17#include "parallelio.h"
18
19class ParallelIoHdf5 : public ParallelIoBase<ParallelIoHdf5> {
21
22 // Public member functions
23 public:
24 // Constructor
25 ParallelIoHdf5(const MString& fileName, MInt fileMode, const MPI_Comm& mpiComm);
26 // Destructor
27 ~ParallelIoHdf5() override;
28
29 // Hdf5-specific methods (private)
30 private:
31 // Static file system-related methods
32 static MBool b_isValidFile(const MString& name, const MPI_Comm& mpiComm);
33 static const MString b_fileExt();
34
35 // File methods
36 void b_open();
37 void b_saveHeader();
40
41
42 // Define mode methods
43 void b_defineArray(maiabd_type type, const MString& name, size_type noDims, size_type* totalCount);
44 void b_defineArray(maiabd_type, const MString&, const MString&, const size_type, const size_type*);
45 void b_defineScalar(maiabd_type type, const MString& name);
46
47 void b_createGroup(const MString& path);
49 size_type b_getNoGroups(const MString& path);
50
51 // Inquiry methods
52 MBool b_hasDataset(const MString& name, const size_type dimension);
53 MBool b_hasDataset(const MString& name, const MString& path);
54 MBool b_hasObject(const MString& path);
55 MInt b_getDatasetType(const MString& name);
56 void b_getDatasetNames(std::vector<MString>& names, const size_type dimension);
57 void b_getDatasetNames(std::vector<MString>& names, const MString& path);
58 void b_getGroupNames(std::vector<MString>& names, const MString& path);
60 size_type b_getDatasetNoDims(const MString& path, const MString& name);
61 size_type b_getDatasetSize(const MString& name, const size_type dimensionId = 0);
62 void b_getDatasetSize(const MString& name, const MString& path, size_type noDims, size_type* data);
63 MBool b_hasAttribute(const MString& name, const MString& path = "");
64 MInt b_getAttributeType(const MString& name, const MString& datasetName = "");
65
66 // Data mode methods
67 template <class T>
68 void b_writeArray(const T* array, const MString& name, const size_type noDims, MPI_Offset* start, MPI_Offset* count,
69 MPI_Offset memoryStride, const size_type noChunks, MPI_Offset diskStride);
70 template <class T>
71 void b_writeArray(const T*, const MString&, const MString&, const size_type, const size_type*, const size_type*,
72 const size_type*);
73 template <class T>
74 void b_writeScalar(T scalar, const MString& name);
75
76 template <class T>
77 void b_readArray(T*, const MString&, const size_type, MPI_Offset*, MPI_Offset*, MPI_Offset, const size_type,
78 MPI_Offset);
79 template <class T>
80 void b_readArray(T*, const MString&, const MString&, const size_type, const size_type*, const size_type*);
81
82 template <class T>
83 void b_readScalar(T* scalar, const MString& name);
84
85 // Attribute methods
86 template <class T>
87 void b_setAttribute(const T* value, const MString& name, const MString& datasetName = "",
88 const size_type totalCount = 1);
89 template <class T>
90 void b_createAttribute(const T* value, const MString& name, const MString& datasetName, hid_t type_id,
91 const size_type totalCount);
92 template <class T>
93 void b_getAttribute(T* value, const MString& name, const MString& datasetName = "", const size_type totalCount = 1);
94 void b_getAttributeCount(const MString& name, size_type* totalCount, const MString& datasetName = "");
95
96 // Auxiliary methods
97 // b_error is static because it is used in b_isValidFile (also static)
98 static void b_error(MInt status, const MString& name, const MString& location);
99 static void b_warning(MInt status, const MString& name, const MString& location);
100
101 // Hdf5-related member variables
102 private:
103 hid_t b_h5Id = -1;
104 // member to manage the I/O operations (collectively)
107};
108#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
MBool b_hasObject(const MString &path)
void b_createAttribute(const T *value, const MString &name, const MString &datasetName, hid_t type_id, const size_type totalCount)
Creates an attribute in the file (generic version).
void b_defineScalar(maiabd_type type, const MString &name)
Defines a scalar in the file. [MPI]
void b_saveHeader()
Adds all additional header information that are needed in the file.
void b_readArray(T *, const MString &, const size_type, MPI_Offset *, MPI_Offset *, MPI_Offset, const size_type, MPI_Offset)
Read array data from file (float version). [MPI]
void b_writeScalar(T scalar, const MString &name)
Writes scalar data to file (Hdf5 version). [MPI]
hid_t b_h5DatasetXferHandle
void b_defineArray(maiabd_type, const MString &, const MString &, const size_type, const size_type *)
void b_writeArray(const 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)
Writes array data to file (Hdf5 version). [MPI]
void b_addAdditionalHeader()
Write additional headers to file (e.g. grid file name, creation date etc.). [MPI]
void b_getAttribute(T *value, const MString &name, const MString &datasetName="", const size_type totalCount=1)
Retrieve an attribute from file (float version).
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.
MBool b_hasAttribute(const MString &name, const MString &path="")
Check if object exists. [MPI]
size_type b_getDatasetNoDims(const MString &name)
Get number of dimensions of a dataset with the given name.
MInt b_getDatasetType(const MString &name)
Returns the data type of a dataset in the file (can be array, multi-D array or scalar).
void b_readScalar(T *scalar, const MString &name)
Read scalar data from file (float version). [MPI]
MInt b_getAttributeType(const MString &name, const MString &datasetName="")
Returns the data type of an attribute in the file.
void b_defineArray(maiabd_type type, const MString &name, size_type noDims, size_type *totalCount)
Defines an array in the file. [MPI]
~ParallelIoHdf5() override
Close open identifiers and release memory. [MPI]
void b_writeAdditionalData()
Write additional data to file. [MPI]
size_type b_getNoDatasets(const MString &path)
Gets the number of datasetes in the given group [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).
static void b_warning(MInt status, const MString &name, const MString &location)
Check the status code of a HDF5 operation and output a meaningful message (Warning)
void b_getGroupNames(std::vector< MString > &names, const MString &path)
Gets the groups in the given group [MPI]
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)....
static const MString b_fileExt()
Returns backend-specific ending of filename (either ".Netcdf" or ".Hdf5")
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_getAttributeCount(const MString &name, size_type *totalCount, const MString &datasetName="")
void b_createGroup(const MString &path)
Creates an HDF5 group [MPI]
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_hasDataset(const MString &name, const size_type dimension)
Check if the file contains a dataset with the given name and number of dimensions.
size_type b_getNoGroups(const MString &path)
Gets the number of groups in the given group [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