MAIA bb96820c
Multiphysics at AIA
|
This class is intended to do all the heavy lifting when it comes to reading and writing "big data" files. It supports multiple backends, i.e. different I/O libraries. The default is "Parallel netCDF". More...
#include <parallelio.h>
Public Types | |
using | size_type = MLong |
Type used for all size- and offset-related values. More... | |
using | maiabd_type = maia::parallel_io::maiabd_type |
Type used for all constants that are defined in maia::parallel_io. More... | |
Public Member Functions | |
MBool | hasDataset (const MString &name, MInt dimension) |
Check if the file contains an dataset with the given name and dimension. More... | |
MBool | hasDataset (const MString &name, const MString &path) |
MBool | hasDataset (const MString &name) |
Check if the file contains an dataset with the given name. More... | |
MBool | hasObject (const MString &path) |
MInt | getDatasetType (const MString &name) |
Returns the data type of an array. More... | |
std::vector< MString > | getDatasetNames (const size_type dimensions=-1) |
Returns a vector with the names of all existing datasets with given dimensionality in the file. More... | |
std::vector< MString > | getDatasetNames (const MString &path) |
std::vector< MString > | getGroupNames (const MString &path) |
size_type | getDatasetNoDims (const MString &name) |
Get the number of dimensions of a dataset with given name. More... | |
size_type | getDatasetNoDims (const MString &name, const MString &path) |
std::vector< size_type > | getArrayDims (const MString &name) |
Get the lengths of all dimensions of a dataset with given name. More... | |
size_type | getArraySize (const MString &name, const size_type dimensionId=0) |
Get the length of an array in the file. More... | |
void | getArraySize (const MString &name, const MString &path, size_type *datasetSize) |
MBool | hasAttribute (const MString &name, const MString &path="") |
Check if a given attribute exists in the file. More... | |
maiabd_type | getAttributeType (const MString &name, const MString &datasetName="") |
Returns the data type of an attribute in the file. More... | |
void | defineArray (maiabd_type type, const MString &name, size_type totalCount) |
Create a new array in the file. More... | |
void | defineArray (maiabd_type type, const MString &name, size_type noDims, size_type *totalCount) |
Create a new array in the file. More... | |
void | defineArray (maiabd_type type, const std::vector< MString > &name, size_type totalCount) |
Create multiple new arrays in the file. More... | |
void | defineArray (maiabd_type type, const MString &datasetName, const MString &name, size_type noDims, size_type *totalCount) |
Creates new array in the file. More... | |
void | defineScalar (maiabd_type type, const MString &name) |
Create a new scalar in the file. More... | |
template<class T > | |
void | writeArray (const T *array, const MString &name, size_type memoryStride=-1, size_type diskStride=-1) |
Write array data to file. [MPI] More... | |
template<class T > | |
void | writeArray (const T *array, const std::vector< MString > &names, size_type memoryStride=-1, size_type diskStride=-1) |
Write array data to file (multi-array version). [MPI] More... | |
template<class T > | |
void | writeArray (const T *array, const MString &path, const MString &name, const size_type noDims, size_type *offset, size_type *localCount) |
Write multidimensional array data to file. [MPI] More... | |
template<class T > | |
void | writeArray (const T *array, const MString &path, const MString &name, const size_type noDims, size_type *offset, size_type *localCount, size_type *ghost) |
Write multidimensional array data to file [ghost cell version]. [MPI] More... | |
template<class T > | |
void | writeScalar (T scalar, const MString &name) |
Write scalar data to file. [MPI] More... | |
template<class T > | |
void | setAttribute (const T &value, const MString &name, const MString &datasetName="") |
Set a file or dataset attribute. [MPI] More... | |
void | setAttribute (const MChar *value, const MString &name, const MString &datasetName="") |
template<class T > | |
void | setAttributes (const T *value, const MString &name, size_type totalCount, const MString &datasetName="") |
template<class T > | |
void | readArray (T *array, const MString &name, size_type memoryStride=-1, size_type diskStride=-1) |
Read array data from file. [MPI] More... | |
template<class T > | |
void | readArray (std::vector< T > &array, const MString &name, size_type memoryStride=-1, size_type diskStride=-1) |
template<class T > | |
void | readArray (T *array, const std::vector< MString > &names, size_type memoryStride=-1, size_type diskStride=-1) |
Read array data from file (multi-array version). [MPI] More... | |
template<class T > | |
void | readArray (T *array, const MString &datasetName, const MString &name, const size_type noDims, const size_type *offset, const size_type *localCount) |
Read array data from file. [MPI] More... | |
template<class T > | |
void | readArray (T *array, const MString &datasetName, const MString &name, const size_type noDims, const size_type *localCount) |
Read array data from file [no offset version]. [MPI] More... | |
template<class T > | |
void | readScalar (T *scalar, const MString &name) |
Read scalar data from file. [MPI] More... | |
template<class T > | |
void | getAttribute (T *value, const MString &name, const MString &datasetName="") |
Retrieve a file or dataset attribute. More... | |
template<class T > | |
void | getAttribute (T *value, const MString &name, const size_type totalCount, const MString &datasetName="") |
void | getAttributeCount (const MString &name, size_type *totalCount, const MString &datasetName="") |
void | setOffset (const size_type localCount, const size_type offset, const size_type noDims, const size_type noChunks) |
Set the local and global counts, as well the local offset for array operations. More... | |
void | setOffset (const size_type localCount, const size_type offset, const size_type noDims) |
Set the local and global counts, as well the local offset for array operations (overload). More... | |
void | setOffset (size_type localCount, size_type offset) |
Set the local and global counts, as well the local offset for array operations (overload). More... | |
Static Public Member Functions | |
static MBool | fileExists (const MString &name, const MPI_Comm &mpiComm) |
Check whether a file exists. [MPI] More... | |
static void | deleteFile (const MString &name) |
Deletes the specified file. [MPI] More... | |
static MBool | isParallelIoFile (const MString &name, const MPI_Comm &mpiComm) |
Check if specified file is a valid ParallelIoHdf5 file. [MPI] More... | |
static MString | fileExt () |
Returns backend-specific ending of filename (either ".Netcdf" or ".Hdf5") More... | |
static void | calcOffset (size_type localCount, size_type *offset, size_type *totalCount, const MPI_Comm &mpiComm, size_type *noChunks=nullptr) |
Calculates the totalCount and offset information needed in setOffset(). [MPI] More... | |
Protected Member Functions | |
ParallelIoBase (MString fileName, MInt fileMode, const MPI_Comm &mpiComm) | |
Creates a new object to read and write big data files. [MPI] More... | |
virtual | ~ParallelIoBase () |
empty destructor. [MPI] More... | |
void | validateType (const maiabd_type type) const |
Auxiliary method that aborts MAIA if type is not a valid type constant for ParallelIo. More... | |
Protected Attributes | |
const MString | m_fileName {} |
const MInt | m_fileMode = -1 |
MBool | m_isHeaderSaved = false |
MInt | m_domainId = -1 |
MInt | m_noDomains = -1 |
MPI_Comm | m_mpiComm {} |
MPI_Info | m_mpiInfo = MPI_INFO_NULL |
MBool | m_offsetsAreSet = false |
std::vector< size_type > | m_localCount = {-1} |
std::vector< size_type > | m_offset = {-1} |
size_type | m_noChunks = -1 |
MFloat | m_creationTime = MFloatNaN |
std::set< MString > | m_unwrittenArrays |
std::set< MString > | m_unwrittenScalars |
std::set< MString > | m_writtenArrays |
std::set< MString > | m_writtenScalars |
std::set< MString > | m_readArrays |
std::set< MString > | m_readScalars |
Static Protected Attributes | |
static const size_type | s_maxChunkSize = 2147483645 |
static const MPI_Datatype | s_mpiSizeType = maia::type_traits<ParallelIoBase<Backend>::size_type>::mpiType() |
Private Member Functions | |
Backend & | derived () |
const Backend & | derived () const |
As of fall 2015 there are two backends available: Parallel netCDF and HDF5.
The main intention of the class is : complete separation of I/O calls from the user, so that nobody has to deal with the different ways of how to use (Parallel) NetCDF, HDF5, MPI I/O etc. The class also offers an easy way to make sure that data files (i.e. restart files, solution files) and grid files always have the correct format.
Note: All methods that may be slow or incur a significant overhead since they use MPI communication (either for I/O or other purposes) are marked with [MPI]. While they might not always have this overhead, it is better to assume they have and to use them only when necessary.
Note: Almost all methods in this class are collective, i.e. they have to be called on all ranks within the given MPI communicator to work properly (even if they themselves are not marked with [MPI]). Some methods may fail if you try to call them from a subset of the ranks, but do not rely on a meaningful error.
Note: Since I/O is often an error-prone area, it is usually advisable to check all I/O operations for possible occurring errors. In order to relieve the user from the burden of doing correct and repeated checking of possible error values/message, ParallelIo adheres to the following rule:
ALL ERRORS ARE FATAL!
Exceptions to this rule are marked as such. This has the advantage that there is no need to do extra checking in the production code of MAIA itself on the one hand, and on the other hand it is assured that it is not possible that I/O operations leave MAIA in an undefined state if an error occurred.
Definition at line 101 of file parallelio.h.
using ParallelIoBase< Backend >::maiabd_type = maia::parallel_io::maiabd_type |
Definition at line 128 of file parallelio.h.
using ParallelIoBase< Backend >::size_type = MLong |
This type will be used internally and externally for all values that can be used as size values or offsets. Examples are 'count', 'size', 'stride' etc. The usage of this type is to ensure that this type may be changed with little to no necessary changes within the code. A change might be necessary e.g. if the previous type's range becomes to small to hold all relevant size information.
Note: If you change this type, you also MUST check if the static member s_mpiSizeType needs to be changed as well, since it is used as the MPI type to communicate values of the 'size_type' kind.
Definition at line 123 of file parallelio.h.
|
protected |
[in] | fileMode | The file mode can be either maia::parallel_io::PIO_CREATE, maia::parallel_io::PIO_APPEND, or maia::parallel_io::PIO_READ. |
[in] | mpiComm | The MPI communicator that should be used to open/create the file. |
Definition at line 405 of file parallelio.h.
|
protectedvirtual |
Definition at line 431 of file parallelio.h.
|
static |
[in] | localCount | Number of elements in an array on the current domain. |
[out] | offset | Starting offset in the file, from which this domain should start writing/reading (is zero if localCount=0). |
[out] | totalCount | Total number of elements in an array across all domains. |
[in] | mpiComm | MPI communicator for which the totalCount and the offset should be calculated. |
[in] | noChunks | Number of Chunks for chunked IO. |
This must be a collective call from all ranks in mpiComm.
Definition at line 1927 of file parallelio.h.
void ParallelIoBase< Backend >::defineArray | ( | maiabd_type | type, |
const MString & | datasetName, | ||
const MString & | name, | ||
size_type | noDims, | ||
size_type * | totalCount | ||
) |
[in] | type | Data type of the arrays (either maia::parallel_io::PIO_FLOAT or maia::parallel_io::PIO_INT). |
[in] | Name | of the array |
[in] | Number | of dimensions of array |
[in] | Size | of array in each space direction |
Definition at line 891 of file parallelio.h.
void ParallelIoBase< Backend >::defineArray | ( | maiabd_type | type, |
const MString & | name, | ||
size_type | noDims, | ||
size_type * | totalCount | ||
) |
[in] | type | Data type of the array (either maia::parallel_io::PIO_FLOAT, maia::parallel_io::PIO_INT, maia::parallel_io::PIO_LONG or maia::parallel_io::PIO_STRING or maia::parallel_io::PIO_UCHAR). |
[in] | name | Name of the array (may not be empty). |
[in] | noDims | Number of array dimensions |
[in] | totalCount | Total size of the array for all dimensions (i.e. combined size of all involved domains). |
This method may not be called after the first call to writeArray() or writeScalar().
Definition at line 812 of file parallelio.h.
void ParallelIoBase< Backend >::defineArray | ( | maiabd_type | type, |
const MString & | name, | ||
size_type | totalCount | ||
) |
[in] | type | Data type of the array (either maia::parallel_io::PIO_FLOAT, maia::parallel_io::PIO_INT, maia::parallel_io::PIO_LONG or maia::parallel_io::PIO_STRING or maia::parallel_io::PIO_UCHAR). |
[in] | name | Name of the array (may not be empty). |
[in] | totalCount | Total size of the array (i.e. combined size of all involved domains). |
This method may not be called after the first call to writeArray() or writeScalar().
Definition at line 783 of file parallelio.h.
void ParallelIoBase< Backend >::defineArray | ( | maiabd_type | type, |
const std::vector< MString > & | names, | ||
size_type | totalCount | ||
) |
[in] | type | Data type of the arrays (either maia::parallel_io::PIO_FLOAT or maia::parallel_io::PIO_INT). |
[in] | names | Names of the arrays (may not be empty). |
[in] | totalCount | Total size of each array (i.e. combined size of all involved domains). |
This method may not be called after the first call to writeArray() or writeScalar().
Definition at line 865 of file parallelio.h.
void ParallelIoBase< Backend >::defineScalar | ( | maiabd_type | type, |
const MString & | name | ||
) |
[in] | type | Data type of the scalar (either maia::parallel_io::PIO_FLOAT or maia::parallel_io::PIO_INT). |
[in] | name | Name of the scalar (may not be empty). |
This method may not be called after the first call to writeArray() or writeScalar().
Definition at line 933 of file parallelio.h.
|
static |
This method does not check if the file exists before attempting to delete it.
Definition at line 346 of file parallelio.h.
|
inlineprivate |
Definition at line 104 of file parallelio.h.
|
inlineprivate |
Definition at line 105 of file parallelio.h.
|
static |
[in] | mpiComm | The MPI communicator to be used to check the file existence. |
This must be a collective call from all ranks in mpiComm.
Definition at line 318 of file parallelio.h.
|
static |
Definition at line 379 of file parallelio.h.
std::vector< typename ParallelIoBase< Backend >::size_type > ParallelIoBase< Backend >::getArrayDims | ( | const MString & | name | ) |
[in] | name | Name of dataset to check |
[out] | dims | Contains the lengths of all dimensions of the dataset |
Definition at line 633 of file parallelio.h.
void ParallelIoBase< Backend >::getArraySize | ( | const MString & | name, |
const MString & | path, | ||
size_type * | datasetSize | ||
) |
Definition at line 688 of file parallelio.h.
ParallelIoBase< Backend >::size_type ParallelIoBase< Backend >::getArraySize | ( | const MString & | name, |
const size_type | dimensionId = 0 |
||
) |
[in] | name | Name of the array to check. |
[in] | dimensionId | Dimension of the array to check (default 0). |
Definition at line 667 of file parallelio.h.
void ParallelIoBase< Backend >::getAttribute | ( | T * | value, |
const MString & | name, | ||
const MString & | datasetName = "" |
||
) |
T | Data type of the attribute (can be either MFloat, MInt, or MString). |
[in] | value | Value of the attribute. |
[in] | name | Name of the attribute. |
[in] | datasetName | Reads the attribute from this dataset (i.e. array or scalar). If empty, retrieve as file attribute. |
Definition at line 1788 of file parallelio.h.
void ParallelIoBase< Backend >::getAttribute | ( | T * | value, |
const MString & | name, | ||
const size_type | totalCount, | ||
const MString & | datasetName = "" |
||
) |
Definition at line 1795 of file parallelio.h.
void ParallelIoBase< Backend >::getAttributeCount | ( | const MString & | name, |
size_type * | totalCount, | ||
const MString & | datasetName = "" |
||
) |
Definition at line 1802 of file parallelio.h.
MInt ParallelIoBase< Backend >::getAttributeType | ( | const MString & | name, |
const MString & | datasetName = "" |
||
) |
[in] | name | Name of the attribute to be checked. |
[in] | datasetName | Name of the dataset for which the attribute is checked. If empty, a file attribute is checked. |
Calls hasAttribute() first to check if the attribute exists and aborts fatally if not.
Definition at line 749 of file parallelio.h.
std::vector< MString > ParallelIoBase< Backend >::getDatasetNames | ( | const MString & | path | ) |
Definition at line 565 of file parallelio.h.
std::vector< MString > ParallelIoBase< Backend >::getDatasetNames | ( | const size_type | dimension = -1 | ) |
[in] | dimension | Return datasets with given dimensionality (default -1) |
[out] | names | Contains all found dataset names with the given dimensionality (if any). |
Definition at line 554 of file parallelio.h.
ParallelIoBase< Backend >::size_type ParallelIoBase< Backend >::getDatasetNoDims | ( | const MString & | name | ) |
[in] | name | Name of dataset to check |
[out] | noDims | Number of dimensions of the dataset |
Definition at line 594 of file parallelio.h.
ParallelIoBase< Backend >::size_type ParallelIoBase< Backend >::getDatasetNoDims | ( | const MString & | name, |
const MString & | path | ||
) |
Definition at line 608 of file parallelio.h.
MInt ParallelIoBase< Backend >::getDatasetType | ( | const MString & | name | ) |
[in] | name | Name of the array to check. |
Calls hasDataset() first to check if the array exists and aborts fatally if not.
Definition at line 530 of file parallelio.h.
std::vector< MString > ParallelIoBase< Backend >::getGroupNames | ( | const MString & | path | ) |
Definition at line 575 of file parallelio.h.
MBool ParallelIoBase< Backend >::hasAttribute | ( | const MString & | name, |
const MString & | path = "" |
||
) |
[in] | name | Name of the attribute to be checked. |
[in] | datasetName | Name of the dataset for which the attribute is checked. If empty, a file attribute is checked. |
Definition at line 714 of file parallelio.h.
MBool ParallelIoBase< Backend >::hasDataset | ( | const MString & | name | ) |
[in] | name | The name of the dataset that should be checked. |
Definition at line 487 of file parallelio.h.
MBool ParallelIoBase< Backend >::hasDataset | ( | const MString & | name, |
const MString & | path | ||
) |
Definition at line 496 of file parallelio.h.
MBool ParallelIoBase< Backend >::hasDataset | ( | const MString & | name, |
MInt | dimension | ||
) |
[in] | name | The name of the dataset that should be checked. |
[in] | dimension | Dimension of dataset to match, -1 will match any dataset, 0 will only match scalars, 1 arrays with one dimensions, ... |
Definition at line 467 of file parallelio.h.
MBool ParallelIoBase< Backend >::hasObject | ( | const MString & | path | ) |
Definition at line 505 of file parallelio.h.
|
static |
[in] | mpiComm | MPI communicator to be used. |
This must be a collective call from all ranks in mpiComm.
Definition at line 364 of file parallelio.h.
|
inline |
Definition at line 195 of file parallelio.h.
void ParallelIoBase< Backend >::readArray | ( | T * | array, |
const MString & | datasetName, | ||
const MString & | name, | ||
const size_type | noDims, | ||
const size_type * | localCount | ||
) |
T | Data type of the array (can be either MFloat or MInt and must match the previous definition of the array). |
[in] | array | Pointer to the data in memory. |
[in] | name | Name of the array in the file. |
[in] | path | Group name in which the array shall be located (works only with HDF5). |
[in] | noDims | Number of dimensions of the array |
[in] | localCount | Size of the array in each dimension, must of of size noDims |
Definition at line 1645 of file parallelio.h.
void ParallelIoBase< Backend >::readArray | ( | T * | array, |
const MString & | datasetName, | ||
const MString & | name, | ||
const size_type | noDims, | ||
const size_type * | offset, | ||
const size_type * | localCount | ||
) |
T | Data type of the array (can be either MFloat or MInt and must match the previous definition of the array). |
[in] | array | Pointer to the data in memory. |
[in] | name | Name of the array in the file. |
[in] | path | Group name in which the array shall be located (works only with HDF5). |
[in] | noDims | Number of dimensions of the array |
[in] | offset | Offset of the of the array in the dataset in each dimension, must be of size noDims |
[in] | localCount | Size of the array in each dimension, must of of size noDims |
Definition at line 1602 of file parallelio.h.
void ParallelIoBase< Backend >::readArray | ( | T * | array, |
const MString & | name, | ||
size_type | memoryStride = -1 , |
||
size_type | diskStride = -1 |
||
) |
T | Data type of the array (can be either MFloat or MInt and must match the previous definition of the array). |
[in] | array | Pointer to the data in memory. |
[in] | name | Name of the array in the file. |
[in] | memoryStride | Offset between two data values in memory (or -1 to set automatically). |
[in] | diskStride | Offset between two data values on disk (or -1 to set automatically). |
The size and file offset have to be set using setOffset() before this method is called. If no stride is given (or -1), the stride will be set to 1.
Note: This method is the exact counterpart to writeArray().
Example:
Suppose you would like to read in the z-coordinates you previously saved to a file in the example of writeArray(maiabd_type, const MString&, size_type). In this case, you just set the offset as before
setOffset(550, 0); // <-- This line on domain0. We have 550 cells and // start at the beginning. setOffset(450, 550); // <-- This line on domain1. We have 450 cells and // leave a gap of 550.
and then call readArray() with the exact same arguments as you used for writeArray():
readArray(&a_coordinates[2], "coordinates_2", 3);
Definition at line 1523 of file parallelio.h.
void ParallelIoBase< Backend >::readArray | ( | T * | array, |
const std::vector< MString > & | names, | ||
size_type | memoryStride = -1 , |
||
size_type | diskStride = -1 |
||
) |
T | Data type of the arrays (can be either MFloat or MInt and must match the previous definition of the array). |
[in] | array | Pointer to the data of the first array in memory. |
[in] | names | Names of the arrays in the file. |
[in] | memoryStride | Offset between two data values in memory (or -1 to set automatically). |
[in] | diskStride | Offset between two data values on disk (or -1 to set automatically). |
The size and file offset have to be set using setOffset() before this method is called. If no stride is given (or -1), the stride will be set to the number of provided names. A call to this method with n elements in 'names' is equivalent to calling readArray() n times, each time with the data pointer 'array' being advanced by 1.
Example:
Suppose you would like to read in the all coordinates you previously saved to a file in the example of writeArray(maiabd_type, const vector<MString>&, size_type). In this case, you just set the offset as before
setOffset(550, 0); // <-- This line on domain0. We have 550 cells and // start at the beginning. setOffset(450, 550); // <-- This line on domain1. We have 450 cells and // leave a gap of 550.
and then call readArray(maiabd_type, const vector<MString>&, size_type) with the exact same arguments as you used for writeArray(maiabd_type, const vector<MString>&, size_type):
readArray(&cells[0].m_coordinates[0], names);
Definition at line 1711 of file parallelio.h.
void ParallelIoBase< Backend >::readScalar | ( | T * | scalar, |
const MString & | name | ||
) |
T | Data type of the scalar (can be either MFloat or MInt and must match the previous definition of the scalar). |
[in] | scalar | Value that should be read from file. |
[in] | name | Name of the scalar in the file. |
Note: This method is the exact counterpart to writeScalar().
Definition at line 1758 of file parallelio.h.
void ParallelIoBase< Backend >::setAttribute | ( | const MChar * | value, |
const MString & | name, | ||
const MString & | datasetName = "" |
||
) |
Definition at line 1445 of file parallelio.h.
void ParallelIoBase< Backend >::setAttribute | ( | const T & | value, |
const MString & | name, | ||
const MString & | datasetName = "" |
||
) |
T | Data type of the attribute (can be either MFloat, MInt, or MString). |
[in] | value | Value of the attribute. |
[in] | name | Name of the attribute (must not be empty). |
[in] | datasetName | Attaches the attribute to this dataset (i.e. array or scalar). If empty, use as file attribute. |
You can add attributes to datasets as soon as they were defined with defineArray()/ defineScalar(). You can also use this to overwrite an attribute of the same name. When you call this method before the first read or write operation, there is no MPI overhead.
Definition at line 1438 of file parallelio.h.
void ParallelIoBase< Backend >::setAttributes | ( | const T * | value, |
const MString & | name, | ||
size_type | totalCount, | ||
const MString & | datasetName = "" |
||
) |
Definition at line 1452 of file parallelio.h.
void ParallelIoBase< Backend >::setOffset | ( | const size_type | localCount, |
const size_type | offset, | ||
const size_type | noDims | ||
) |
Definition at line 1889 of file parallelio.h.
void ParallelIoBase< Backend >::setOffset | ( | const size_type | localCount, |
const size_type | offset, | ||
const size_type | noDims, | ||
const size_type | noChunks | ||
) |
[in] | localCount | Number of elements in an array on the current domain. |
[in] | offset | Starting offset in the file, from which this domain should start writing/reading. |
[in] | noDims | number of dimensions of the array for more than one dimension localCount and offset refer to the first dimension, for all other dimensions the count will be set to the total dimension length, multidimensional counts and offsets are not yet supported |
[in] | noChunks | Number of Chunks for chunked IO. |
You have to call this method at least once before your first call to readArray()/writeArray(). The reason for storing the offset and count information with the object is that there are usually only a very limited number of different array sizes in each file, and thus it is more convenient to set the offsets once and write/read multiple times afterwards.
The offset may be calculated by hand (or known beforehand), or may be conveniently calculated using calcOffset().
Definition at line 1840 of file parallelio.h.
void ParallelIoBase< Backend >::setOffset | ( | size_type | localCount, |
size_type | offset | ||
) |
Definition at line 1901 of file parallelio.h.
|
protected |
[in] | type | The type constant (integer) to check for validity. |
Definition at line 1976 of file parallelio.h.
void ParallelIoBase< Backend >::writeArray | ( | const T * | array, |
const MString & | name, | ||
size_type | memoryStride = -1 , |
||
size_type | diskStride = -1 |
||
) |
T | Data type of the array (can be either MFloat or MInt and must match the previous definition of the array). |
[in] | array | Pointer to the data in memory. |
[in] | name | Name of the array in the file. |
[in] | memoryStride | Offset between two data values in memory (or -1 to set automatically). |
[in] | diskStride | Offset between two data values on disk (or -1 to set automatically). |
The size and file offset have to be set using setOffset() before this method is called. If no stride is given (or -1), the stride will be set to 1.
Example:
Suppose you would like to store all z-coordinates of your cells in a file, and assign them the name 'coordinates_2'. You have a total of 1000 cells which are split across two domains, domain0 with 550 cells and domain1 with the remaing 450.
Then you first have to define the array with a call
defineArray(PIO_FLOAT, "coordinates_2", 1000);
Note that '1000' denoting the total count, i.e. the number of values you would like to write from all domains (or, alternatively worded, this is the size the dataset will have in the file). If you do not know the total count yet, you can use calcOffset() to conveniently calculate it.
Next you have to set the offset, so each processor knows where it should start writing. In this case its rather simple: domain0 should write its first value at the beginning of the file, and domain1 should leave a gap of 550 values for domain0. Thus you have to call setOffset() once on each domain, i.e.
setOffset(550, 0); // <-- This line on domain0. We have 550 cells and // start at the beginning. setOffset(450, 550); // <-- This line on domain1. We have 450 cells and // leave a gap of 550.
Now you can start writing. For this you need to know how your data is distributed in memory. In MAIA, most cell data is stored contiguously in memory, grouped by cell. This means that for 3D, the coordinates are layed out in memory as (noted as cellId:direction)
| 0:x | 0:y | 0:z | 1:x | 1:y | 1:z | 2:x | 2:y | 2:z | 3:x | 3:y | ...
For the writeArray() methods you always need to get a pointer to the first value in memory that should be written. The memory address of the z-coordinate of the first cell can be obtained using '&cells[0].m_coordinates[2]'. Therefore, in order to write your data you make the same call to writeArray() on both domains:
writeArray(&a_coordinates[2], "coordinates_2", 3);
What is the '3' doing there? Well, since the data values in memory are not contiguous (we only want to write the z-coordinates, i.e. we only need every third value), you need to provide a stride. The stride always refers to the layout in memory - in the file the values will always be contiguous.
Finally, when we check the file, we will see the following data in 'coordinates_2':
| 0:z | 1:z | 2:z | 3:z | 4:z | ... | 998:z | 999:z |
Now, that wasn't so hard, was it?
Definition at line 1046 of file parallelio.h.
void ParallelIoBase< Backend >::writeArray | ( | const T * | array, |
const MString & | path, | ||
const MString & | name, | ||
const size_type | noDims, | ||
size_type * | offset, | ||
size_type * | localCount | ||
) |
T | Data type of the array (can be either MFloat or MInt and must match the previous definition of the array). |
[in] | array | Pointer to the data in memory. |
[in] | path | Group name in which the array shall be located (works only with HDF5). |
[in] | name | Name of the array in the file. |
[in] | noDims | Number of dimensions of the array |
[in] | offset | Offset of the of the array in the dataset in each dimension, must be of size noDims |
[in] | localCount | Size of the array in each dimension, must of of size noDims |
Definition at line 1255 of file parallelio.h.
void ParallelIoBase< Backend >::writeArray | ( | const T * | array, |
const MString & | path, | ||
const MString & | name, | ||
const size_type | noDims, | ||
size_type * | offset, | ||
size_type * | localCount, | ||
size_type * | ghost | ||
) |
T | Data type of the array (can be either MFloat or MInt and must match the previous definition of the array). |
[in] | array | Pointer to the data in memory. |
[in] | path | Group name in which the array shall be located (works only with HDF5). |
[in] | name | Name of the array in the file. |
[in] | noDims | Number of dimensions of the array |
[in] | offset | Offset of the of the array in the dataset in each dimension, must be of size noDims |
[in] | localCount | Size of the array in each dimension, must of of size noDims |
[in] | ghost | Number of ghost cells in memory preceding the actual values |
This method is specifically for storing two- or three-dimensional arrays from the structured solver. Typically the in-memory array has a number of ghost cells at the beginning in every dimension (and at the end) which should be added to the offset when reading from memory, but should be ignored for the dataset in the file that you're writing to. That is, the ghost values are used for the access of the data in memory but ignored when they are being written to the actual dataset in the file.
Definition at line 1292 of file parallelio.h.
void ParallelIoBase< Backend >::writeArray | ( | const T * | array, |
const std::vector< MString > & | names, | ||
size_type | memoryStride = -1 , |
||
size_type | diskStride = -1 |
||
) |
T | Data type of the arrays (can be either MFloat or MInt and must match the previous definition of the array). |
[in] | array | Pointer to the data of the first array in memory. |
[in] | names | Names of the arrays in the file. |
[in] | memoryStride | Offset between two data values in memory (or -1 to set automatically). |
[in] | diskStride | Offset between two data values on disk (or -1 to set automatically). |
The size and file offset have to be set using setOffset() before this method is called. If no stride is given (or -1), the stride will be set to the number of provided names. A call to this method with n elements in 'names' is equivalent to calling writeArray() n times, each time with the data pointer 'array' being advanced by 1.
Example:
Suppose you would like to write all coordinates of all cells to a file, and store them as separate variables 'coordinates_0' ... 'coordinates_2' (in 3D). Similarly to the example in writeArray(), you have two domains with 550 cells on domain0, and 450 cells on domain1.
First you have to call defineArray() and setOffset(). You call both with the same arguments as described in writeArray(). The only difference is that instead of defining each array individually, you can call defineArray(maiabd_type, const vector<MString>&, size_type) with a vector that holds your names (in our case the aforementioned 'coordinates_0', 'coordinates_1', and 'coordinates_2'), which we call 'names' in this example:
defineArray(PIO_FLOAT, names, 1000);
This is fully equivalent to writing
defineArray(PIO_FLOAT, 'coordinates_0', 1000); defineArray(PIO_FLOAT, 'coordinates_1', 1000); defineArray(PIO_FLOAT, 'coordinates_2', 1000);
Next, to set the offset, you call setOffset():
setOffset(550, 0); // <-- This line on domain0. We have 550 cells and // start at the beginning. setOffset(450, 550); // <-- This line on domain1. We have 450 cells and // leave a gap of 550.
Finally, we can write the data. Using the vector-of-strings version of writeArray(maiabd_type, const vector<MString>&, size_type), we only need one call:
writeArray(&a_coordinates[0], names, 3);
This is fully equivalent to writing
writeArray(&a_coordinates[0], 'coordinates_0', 3); writeArray(&a_coordinates[1], 'coordinates_1', 3); writeArray(&a_coordinates[2], 'coordinates_2', 3);
but more convenient :-). It gets even better: If you have want to write all components of a multi-component array from memory to the file (as in this example), writeArray(maiabd_type, const vector<MString>&, size_type) can do the stride calculation automatically for you. So we could just leave it out and write
writeArray(&cells[0].m_coordinates[0], names);
instead.
Definition at line 1204 of file parallelio.h.
void ParallelIoBase< Backend >::writeScalar | ( | T | scalar, |
const MString & | name | ||
) |
T | Data type of the scalar (can be either MFloat or MInt and must match the previous definition of the scalar). |
[in] | scalar | Value that should be written to file. |
[in] | name | Name of the scalar in the file. |
Definition at line 1354 of file parallelio.h.
|
protected |
Definition at line 250 of file parallelio.h.
|
protected |
Definition at line 236 of file parallelio.h.
|
protected |
Definition at line 231 of file parallelio.h.
|
protected |
Definition at line 230 of file parallelio.h.
|
protected |
Definition at line 232 of file parallelio.h.
|
protected |
Definition at line 245 of file parallelio.h.
|
protected |
Definition at line 240 of file parallelio.h.
|
protected |
Definition at line 241 of file parallelio.h.
|
protected |
Definition at line 247 of file parallelio.h.
|
protected |
Definition at line 237 of file parallelio.h.
|
protected |
Definition at line 246 of file parallelio.h.
|
protected |
Definition at line 244 of file parallelio.h.
|
protected |
Definition at line 258 of file parallelio.h.
|
protected |
Definition at line 259 of file parallelio.h.
|
protected |
Definition at line 254 of file parallelio.h.
|
protected |
Definition at line 255 of file parallelio.h.
|
protected |
Definition at line 256 of file parallelio.h.
|
protected |
Definition at line 257 of file parallelio.h.
|
staticprotected |
Definition at line 265 of file parallelio.h.
|
staticprotected |
Definition at line 270 of file parallelio.h.