MAIA bb96820c
Multiphysics at AIA
Loading...
Searching...
No Matches
ParallelIoBase< Backend > Class Template Reference

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>

Inheritance diagram for ParallelIoBase< Backend >:
[legend]
Collaboration diagram for ParallelIoBase< Backend >:
[legend]

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< MStringgetDatasetNames (const size_type dimensions=-1)
 Returns a vector with the names of all existing datasets with given dimensionality in the file. More...
 
std::vector< MStringgetDatasetNames (const MString &path)
 
std::vector< MStringgetGroupNames (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_typegetArrayDims (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_typem_localCount = {-1}
 
std::vector< size_typem_offset = {-1}
 
size_type m_noChunks = -1
 
MFloat m_creationTime = MFloatNaN
 
std::set< MStringm_unwrittenArrays
 
std::set< MStringm_unwrittenScalars
 
std::set< MStringm_writtenArrays
 
std::set< MStringm_writtenScalars
 
std::set< MStringm_readArrays
 
std::set< MStringm_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
 

Detailed Description

template<class Backend>
class ParallelIoBase< Backend >
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-03-27

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.

Member Typedef Documentation

◆ maiabd_type

template<class Backend >
using ParallelIoBase< Backend >::maiabd_type = maia::parallel_io::maiabd_type

Definition at line 128 of file parallelio.h.

◆ size_type

template<class Backend >
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.

Constructor & Destructor Documentation

◆ ParallelIoBase()

template<class Backend >
ParallelIoBase< Backend >::ParallelIoBase ( MString  fileName,
MInt  fileMode,
const MPI_Comm &  mpiComm 
)
protected
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2013-03-26
Author
(Modified by: )Ramandeep Jain (HiWi) raman.nosp@m.deep.nosp@m.jain@.nosp@m.gmai.nosp@m.l.com, Konstantin Froehlich
Date
2015-03-01
Parameters
[in]fileModeThe file mode can be either maia::parallel_io::PIO_CREATE, maia::parallel_io::PIO_APPEND, or maia::parallel_io::PIO_READ.
[in]mpiCommThe MPI communicator that should be used to open/create the file.

Definition at line 405 of file parallelio.h.

406 : m_fileName(std::move(fileName)), m_fileMode(fileMode), m_mpiComm(mpiComm) {
407 TRACE();
408
409 using namespace maia::parallel_io;
410 // Make sure that the file mode is supported
411 if(fileMode != PIO_CREATE && fileMode != PIO_APPEND && fileMode != PIO_REPLACE && fileMode != PIO_READ) {
412 TERMM(1, "File mode must be one of PIO_CREATE/PIO_APPEND/PIO_REPLACE/PIO_READ!");
413 }
414
415 // Set domain id and number of domains
416 MPI_Comm_rank(m_mpiComm, &m_domainId);
417 MPI_Comm_size(m_mpiComm, &m_noDomains);
418
419 // Store creation time stamp
421}
const MString m_fileName
Definition: parallelio.h:230
MPI_Comm m_mpiComm
Definition: parallelio.h:240
MFloat m_creationTime
Definition: parallelio.h:250
const MInt m_fileMode
Definition: parallelio.h:231
MFloat wallTime()
Definition: functions.h:80

◆ ~ParallelIoBase()

template<class Backend >
ParallelIoBase< Backend >::~ParallelIoBase
protectedvirtual
Author
Konstantin Froehlich
Date
2013-03-26

Definition at line 431 of file parallelio.h.

431 {
432 TRACE();
433
434 // Log the total file lifetime
435 const MFloat lifetime = wallTime() - m_creationTime;
436
437 const MInt maxLineLength = 256;
438 MChar b[maxLineLength];
439 snprintf(b, maxLineLength, "=== PARALLELIO FILE LIFETIME: %-35s | %.4e s |\n", m_fileName.c_str(), lifetime);
440 if(m_domainId == 0) {
441 std::cerr << b;
442 }
443 m_log << b;
444}
InfoOutFile m_log
int32_t MInt
Definition: maiatypes.h:62
double MFloat
Definition: maiatypes.h:52
char MChar
Definition: maiatypes.h:56

Member Function Documentation

◆ calcOffset()

template<class Backend >
void ParallelIoBase< Backend >::calcOffset ( size_type  localCount,
size_type offset,
size_type totalCount,
const MPI_Comm &  mpiComm,
size_type noChunks = nullptr 
)
static
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-09
Parameters
[in]localCountNumber of elements in an array on the current domain.
[out]offsetStarting offset in the file, from which this domain should start writing/reading (is zero if localCount=0).
[out]totalCountTotal number of elements in an array across all domains.
[in]mpiCommMPI communicator for which the totalCount and the offset should be calculated.
[in]noChunksNumber of Chunks for chunked IO.

This must be a collective call from all ranks in mpiComm.

Definition at line 1927 of file parallelio.h.

1928 {
1929 TRACE();
1930
1931 // Ensure reasonable parameter values
1932 if(localCount < 0) {
1933 TERMM(1, "The local size may not be less than zero (was:" + std::to_string(localCount) + ")!");
1934 }
1935
1936 // Get rank and size
1937 MInt domainId, noDomains;
1938 MPI_Comm_rank(mpiComm, &domainId);
1939 MPI_Comm_size(mpiComm, &noDomains);
1940
1941 // Calculate offset and total size values
1942 MPI_Exscan(&localCount, offset, 1, s_mpiSizeType, MPI_SUM, mpiComm, AT_, "localCount", "offset");
1943 if(domainId == 0) {
1944 offset[0] = 0;
1945 }
1946
1947 totalCount[0] = offset[0] + localCount;
1948 MPI_Bcast(totalCount, 1, s_mpiSizeType, (noDomains - 1), mpiComm, AT_, "totalCount");
1949
1950 // If a domain has zero size set the offset to zero such that reading/writing
1951 // can not result in an error caused by exceeding the array size
1952 if(localCount == 0) {
1953 offset[0] = 0;
1954 }
1955
1956 // If noChunks is not a null pointer, calculate & set number of chunks
1957 if(noChunks != nullptr) {
1958 // Divide total count by maximum chunk size
1959 noChunks[0] = totalCount[0] / s_maxChunkSize;
1960
1961 // If there is a remainder for the division, increment number by one
1962 if(totalCount[0] % s_maxChunkSize > 0) {
1963 noChunks[0] += 1;
1964 }
1965 }
1966}
static const size_type s_maxChunkSize
Definition: parallelio.h:265
static const MPI_Datatype s_mpiSizeType
Definition: parallelio.h:270
int MPI_Exscan(const void *sendbuf, void *recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm, const MString &name, const MString &sndvarname, const MString &rcvvarname)
same as MPI_Exscan
int MPI_Bcast(void *buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm, const MString &name, const MString &varname)
same as MPI_Bcast

◆ defineArray() [1/4]

template<class Backend >
void ParallelIoBase< Backend >::defineArray ( maiabd_type  type,
const MString datasetName,
const MString name,
size_type  noDims,
size_type totalCount 
)
Author
Marian Albers <m.alb.nosp@m.ers@.nosp@m.aia.r.nosp@m.wth-.nosp@m.aache.nosp@m.n.de
Date
2022-09-25
Parameters
[in]typeData type of the arrays (either maia::parallel_io::PIO_FLOAT or maia::parallel_io::PIO_INT).
[in]Nameof the array
[in]Numberof dimensions of array
[in]Sizeof array in each space direction

Definition at line 891 of file parallelio.h.

892 {
893 TRACE();
894
895 using namespace maia::parallel_io;
896
897 // Prevent adding data in read-only mode
898 if(m_fileMode == PIO_READ) {
899 mTerm(1, AT_, "Cannot add data in read-only mode!");
900 }
901
902#ifdef MAIA_EXTRA_DEBUG
903 m_unwrittenArrays.insert(name);
904#endif
905
906 // Ensure reasonable parameter values
907 validateType(type);
908 if(name.empty()) {
909 TERMM(1, "Invalid name! Name must not be empty.");
910 }
911 if(noDims <= 0) {
912 TERMM(1, "Invalid number of dimensions! Must at least be 1.");
913 }
914
915 derived().b_defineArray(type, datasetName, name, noDims, totalCount);
916}
Backend & derived()
Definition: parallelio.h:104
void validateType(const maiabd_type type) const
Auxiliary method that aborts MAIA if type is not a valid type constant for ParallelIo.
Definition: parallelio.h:1976
std::set< MString > m_unwrittenArrays
Definition: parallelio.h:254
void mTerm(const MInt errorCode, const MString &location, const MString &message)
Definition: functions.cpp:29

◆ defineArray() [2/4]

template<class Backend >
void ParallelIoBase< Backend >::defineArray ( maiabd_type  type,
const MString name,
size_type  noDims,
size_type totalCount 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-08
Parameters
[in]typeData 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]nameName of the array (may not be empty).
[in]noDimsNumber of array dimensions
[in]totalCountTotal 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.

813 {
814 TRACE();
815
816#ifdef DISABLE_OUTPUT
817 return;
818#endif
819
820 using namespace maia::parallel_io;
821
822 // Prevent adding data in read-only mode
823 if(m_fileMode == PIO_READ) {
824 mTerm(1, AT_, "Cannot add data in read-only mode!");
825 }
826
827 // Abort if header was already written
828 if(m_isHeaderSaved) {
829 mTerm(1, AT_, "Cannot add new data after header was written!");
830 }
831
832#ifdef MAIA_EXTRA_DEBUG
833 m_unwrittenArrays.insert(name);
834#endif
835
836 // Ensure reasonable parameter values
837 validateType(type);
838 if(name.empty()) {
839 TERMM(1, "Invalid name! Name must not be empty.");
840 }
841 if(noDims <= 0) {
842 TERMM(1, "Invalid number of dimensions! Must at least be 1.");
843 }
844
845 derived().b_defineArray(type, name, noDims, totalCount);
846}
MBool m_isHeaderSaved
Definition: parallelio.h:232

◆ defineArray() [3/4]

template<class Backend >
void ParallelIoBase< Backend >::defineArray ( maiabd_type  type,
const MString name,
size_type  totalCount 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-08
Parameters
[in]typeData 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]nameName of the array (may not be empty).
[in]totalCountTotal 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.

783 {
784 TRACE();
785
786#ifdef DISABLE_OUTPUT
787 return;
788#endif
789
790 defineArray(type, name, 1, &totalCount);
791}
void defineArray(maiabd_type type, const MString &name, size_type totalCount)
Create a new array in the file.
Definition: parallelio.h:783

◆ defineArray() [4/4]

template<class Backend >
void ParallelIoBase< Backend >::defineArray ( maiabd_type  type,
const std::vector< MString > &  names,
size_type  totalCount 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-08
Parameters
[in]typeData type of the arrays (either maia::parallel_io::PIO_FLOAT or maia::parallel_io::PIO_INT).
[in]namesNames of the arrays (may not be empty).
[in]totalCountTotal 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.

865 {
866 TRACE();
867
868#ifdef DISABLE_OUTPUT
869 return;
870#endif
871
872 for(const auto& name : names) {
873 defineArray(type, name, totalCount);
874 }
875}

◆ defineScalar()

template<class Backend >
void ParallelIoBase< Backend >::defineScalar ( maiabd_type  type,
const MString name 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-08
Parameters
[in]typeData type of the scalar (either maia::parallel_io::PIO_FLOAT or maia::parallel_io::PIO_INT).
[in]nameName 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.

933 {
934 TRACE();
935
936#ifdef DISABLE_OUTPUT
937 return;
938#endif
939
940 using namespace maia::parallel_io;
941
942 // Prevent adding data in read-only mode
943 if(m_fileMode == PIO_READ) {
944 TERMM(1, "Cannot add data in read-only mode!");
945 }
946
947 // Abort if header was already written
948 if(m_isHeaderSaved) {
949 TERMM(1, "Cannot add new data after header was written!");
950 }
951
952#ifdef MAIA_EXTRA_DEBUG
953 m_unwrittenScalars.insert(name);
954#endif
955
956 // Ensure reasonable parameter values
957 validateType(type);
958 if(name.empty()) {
959 TERMM(1, "Invalid name! Name must not be empty.");
960 }
961
962 derived().b_defineScalar(type, name);
963}
std::set< MString > m_unwrittenScalars
Definition: parallelio.h:255

◆ deleteFile()

template<typename Backend >
void ParallelIoBase< Backend >::deleteFile ( const MString name)
static
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2013-04-01

This method does not check if the file exists before attempting to delete it.

Definition at line 346 of file parallelio.h.

346 {
347 TRACE();
348
349 MPI_File_delete(const_cast<MChar*>(name.c_str()), MPI_INFO_NULL);
350}

◆ derived() [1/2]

template<class Backend >
Backend & ParallelIoBase< Backend >::derived ( )
inlineprivate

Definition at line 104 of file parallelio.h.

104{ return static_cast<Backend&>(*this); }

◆ derived() [2/2]

template<class Backend >
const Backend & ParallelIoBase< Backend >::derived ( ) const
inlineprivate

Definition at line 105 of file parallelio.h.

105{ return static_cast<const Backend&>(*this); }

◆ fileExists()

template<typename Backend >
MBool ParallelIoBase< Backend >::fileExists ( const MString name,
const MPI_Comm &  mpiComm 
)
static
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-01
Parameters
[in]mpiCommThe MPI communicator to be used to check the file existence.
Returns
True, if the specified file exists.

This must be a collective call from all ranks in mpiComm.

Definition at line 318 of file parallelio.h.

318 {
319 TRACE();
320
321 // Try to stat file on rank 0 - will return something other than zero of file
322 // does not exist - and broadcast to others
323 int rank, status;
324 MPI_Comm_rank(mpiComm, &rank);
325 if(rank == 0) {
326 struct stat buffer {};
327 status = static_cast<int>(stat(name.c_str(), &buffer) == 0);
328 }
329 MPI_Bcast(&status, 1, MPI_INT, 0, mpiComm, AT_, "status");
330
331 return status != 0;
332}

◆ fileExt()

template<class Backend >
MString ParallelIoBase< Backend >::fileExt
static
Author
Konstantin Froehlich
Date
2015-10-21

Definition at line 379 of file parallelio.h.

379 {
380 TRACE();
381
382 return Backend::b_fileExt();
383}

◆ getArrayDims()

template<class Backend >
std::vector< typename ParallelIoBase< Backend >::size_type > ParallelIoBase< Backend >::getArrayDims ( const MString name)
Author
Ansgar Niemoeller, Konstantin Froehlich
Date
2015-08-19
Parameters
[in]nameName of dataset to check
[out]dimsContains the lengths of all dimensions of the dataset

Definition at line 633 of file parallelio.h.

633 {
634 TRACE();
635
636 // Abort if array does not exist
637 if(!hasDataset(name) || hasDataset(name, 0)) {
638 TERMM(1, "The specified array '" + name + "' does not exist.");
639 }
640
641 // Get number of dimensions
642 size_type noDims = getDatasetNoDims(name);
643
644 std::vector<size_type> dims;
645
646 // Get lengths of all dimensions
647 for(size_type dimId = 0; dimId < noDims; dimId++) {
648 dims.push_back(getArraySize(name, dimId));
649 }
650
651 return dims;
652}
MBool hasDataset(const MString &name, MInt dimension)
Check if the file contains an dataset with the given name and dimension.
Definition: parallelio.h:467
MLong size_type
Type used for all size- and offset-related values.
Definition: parallelio.h:123
size_type getArraySize(const MString &name, const size_type dimensionId=0)
Get the length of an array in the file.
Definition: parallelio.h:667
size_type getDatasetNoDims(const MString &name)
Get the number of dimensions of a dataset with given name.
Definition: parallelio.h:594

◆ getArraySize() [1/2]

template<class Backend >
void ParallelIoBase< Backend >::getArraySize ( const MString name,
const MString path,
size_type datasetSize 
)

Definition at line 688 of file parallelio.h.

688 {
689 TRACE();
690
691 // Abort if dataset does not exist or if its a scalar
692 if(!hasDataset(name, path)) {
693 TERMM(1, "The specified array '" + name + "' does not exist.");
694 }
695
696 size_type noDims = derived().b_getDatasetNoDims(name, path);
697 derived().b_getDatasetSize(name, path, noDims, datasetSize);
698}

◆ getArraySize() [2/2]

template<class Backend >
ParallelIoBase< Backend >::size_type ParallelIoBase< Backend >::getArraySize ( const MString name,
const size_type  dimensionId = 0 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-01
Parameters
[in]nameName of the array to check.
[in]dimensionIdDimension of the array to check (default 0).
Returns
Total number of elements of the array in the given dimension.

Definition at line 667 of file parallelio.h.

668 {
669 TRACE();
670
671 // Abort if dataset does not exist or if its a scalar
672 if(!hasDataset(name) || hasDataset(name, 0)) {
673 TERMM(1, "The specified array '" + name + "' does not exist.");
674 }
675
676 // Abort if dimensionId exceeds number of dimensions
677 size_type noDims = getDatasetNoDims(name);
678 if(dimensionId >= noDims) {
679 TERMM(1, "The specified dimension for array '" + name + "'does not exist.");
680 }
681
682 size_type returnValue = derived().b_getDatasetSize(name, dimensionId);
683
684 return returnValue;
685}

◆ getAttribute() [1/2]

template<class Backend >
template<class T >
void ParallelIoBase< Backend >::getAttribute ( T *  value,
const MString name,
const MString datasetName = "" 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-09
Template Parameters
TData type of the attribute (can be either MFloat, MInt, or MString).
Parameters
[in]valueValue of the attribute.
[in]nameName of the attribute.
[in]datasetNameReads the attribute from this dataset (i.e. array or scalar). If empty, retrieve as file attribute.

Definition at line 1788 of file parallelio.h.

1788 {
1789 TRACE();
1790 getAttribute(value, name, 1, datasetName);
1791}
void getAttribute(T *value, const MString &name, const MString &datasetName="")
Retrieve a file or dataset attribute.
Definition: parallelio.h:1788

◆ getAttribute() [2/2]

template<class Backend >
template<class T >
void ParallelIoBase< Backend >::getAttribute ( T *  value,
const MString name,
const size_type  totalCount,
const MString datasetName = "" 
)

Definition at line 1795 of file parallelio.h.

1796 {
1797 TRACE();
1798 derived().b_getAttribute(value, name, datasetName, totalCount);
1799}

◆ getAttributeCount()

template<class Backend >
void ParallelIoBase< Backend >::getAttributeCount ( const MString name,
size_type totalCount,
const MString datasetName = "" 
)

Definition at line 1802 of file parallelio.h.

1803 {
1804 TRACE();
1805 derived().b_getAttributeCount(name, totalCount, datasetName);
1806}

◆ getAttributeType()

template<class Backend >
MInt ParallelIoBase< Backend >::getAttributeType ( const MString name,
const MString datasetName = "" 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-01
Parameters
[in]nameName of the attribute to be checked.
[in]datasetNameName of the dataset for which the attribute is checked. If empty, a file attribute is checked.
Returns
The data type as defined in the namespace maia::parallel_io.

Calls hasAttribute() first to check if the attribute exists and aborts fatally if not.

Definition at line 749 of file parallelio.h.

749 {
750 TRACE();
751
752 // Abort if attribute does not exist
753 if(!hasAttribute(name, datasetName)) {
754 mTerm(1, AT_, "The specified attribute does not exist.");
755 }
756
757 return derived().b_getAttributeType(name, datasetName);
758}
MBool hasAttribute(const MString &name, const MString &path="")
Check if a given attribute exists in the file.
Definition: parallelio.h:714

◆ getDatasetNames() [1/2]

template<class Backend >
std::vector< MString > ParallelIoBase< Backend >::getDatasetNames ( const MString path)

Definition at line 565 of file parallelio.h.

565 {
566 TRACE();
567
568 std::vector<MString> names;
569 derived().b_getDatasetNames(names, path);
570
571 return names;
572}

◆ getDatasetNames() [2/2]

template<class Backend >
std::vector< MString > ParallelIoBase< Backend >::getDatasetNames ( const size_type  dimension = -1)
Author
Ansgar Niemoeller, , Konstantin Froehlich
Date
2015-08-19
Parameters
[in]dimensionReturn datasets with given dimensionality (default -1)
[out]namesContains all found dataset names with the given dimensionality (if any).

Definition at line 554 of file parallelio.h.

554 {
555 TRACE();
556
557 std::vector<MString> names;
558 derived().b_getDatasetNames(names, dimension);
559
560 return names;
561}

◆ getDatasetNoDims() [1/2]

template<class Backend >
ParallelIoBase< Backend >::size_type ParallelIoBase< Backend >::getDatasetNoDims ( const MString name)
Author
Ansgar Niemoeller, Konstantin Froehlich
Date
2015-08-19
Parameters
[in]nameName of dataset to check
[out]noDimsNumber of dimensions of the dataset

Definition at line 594 of file parallelio.h.

594 {
595 TRACE();
596
597 // Abort if dataset does not exist
598 if(!hasDataset(name)) {
599 TERMM(1, "The specified dataset '" + name + "' does not exist.");
600 }
601
602 size_type noDims = derived().b_getDatasetNoDims(name);
603
604 return noDims;
605}

◆ getDatasetNoDims() [2/2]

template<class Backend >
ParallelIoBase< Backend >::size_type ParallelIoBase< Backend >::getDatasetNoDims ( const MString name,
const MString path 
)

Definition at line 608 of file parallelio.h.

609 {
610 TRACE();
611
612 // Abort if dataset does not exist
613 if(!hasDataset(name, path)) {
614 TERMM(1, "The specified dataset '" + name + "' does not exist.");
615 }
616
617 size_type noDims = derived().b_getDatasetNoDims(name, path);
618
619 return noDims;
620}

◆ getDatasetType()

template<class Backend >
MInt ParallelIoBase< Backend >::getDatasetType ( const MString name)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2013-04-01
Author
Modified by: Ramandeep Jain (HiWi) raman.nosp@m.deep.nosp@m.jain@.nosp@m.gmai.nosp@m.l.com, Konstantin Froehlich
Date
2015-Mar
Parameters
[in]nameName of the array to check.
Returns
The data type as defined in the namespace maia::parallel_io.

Calls hasDataset() first to check if the array exists and aborts fatally if not.

Definition at line 530 of file parallelio.h.

530 {
531 TRACE();
532
533 // Abort if dataset does not exist
534 if(!hasDataset(name)) {
535 TERMM(1, "The specified dataset '" + name + "' does not exist.");
536 }
537
538 return derived().b_getDatasetType(name);
539}

◆ getGroupNames()

template<class Backend >
std::vector< MString > ParallelIoBase< Backend >::getGroupNames ( const MString path)

Definition at line 575 of file parallelio.h.

575 {
576 TRACE();
577
578 std::vector<MString> names;
579 derived().b_getGroupNames(names, path);
580
581 return names;
582}

◆ hasAttribute()

template<class Backend >
MBool ParallelIoBase< Backend >::hasAttribute ( const MString name,
const MString path = "" 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-01
Parameters
[in]nameName of the attribute to be checked.
[in]datasetNameName of the dataset for which the attribute is checked. If empty, a file attribute is checked.
Returns
True, if the attribute exists.

Definition at line 714 of file parallelio.h.

714 {
715 TRACE();
716
717 // If the dataset name is empty, this means a file attribute should be
718 // checked.
719 // Otherwise a dataset attribute (scalar or array) is to be checked - abort if
720 // no dataset can be found.
721 if(!path.empty()) {
722 if(!hasObject(path)) {
723 return false;
724 }
725 }
726
727 MBool returnValue = derived().b_hasAttribute(name, path);
728
729 return returnValue;
730}
MBool hasObject(const MString &path)
Definition: parallelio.h:505
bool MBool
Definition: maiatypes.h:58

◆ hasDataset() [1/3]

template<class Backend >
MBool ParallelIoBase< Backend >::hasDataset ( const MString name)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-01
Parameters
[in]nameThe name of the dataset that should be checked.
Returns
True if an dataset of the given name exists.

Definition at line 487 of file parallelio.h.

487 {
488 TRACE();
489
490 MBool returnValue = derived().b_hasDataset(name, -1);
491
492 return returnValue;
493}

◆ hasDataset() [2/3]

template<class Backend >
MBool ParallelIoBase< Backend >::hasDataset ( const MString name,
const MString path 
)

Definition at line 496 of file parallelio.h.

496 {
497 TRACE();
498
499 MBool returnValue = derived().b_hasDataset(name, path);
500
501 return returnValue;
502}

◆ hasDataset() [3/3]

template<class Backend >
MBool ParallelIoBase< Backend >::hasDataset ( const MString name,
MInt  dimension 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-01
Parameters
[in]nameThe name of the dataset that should be checked.
[in]dimensionDimension of dataset to match, -1 will match any dataset, 0 will only match scalars, 1 arrays with one dimensions, ...
Returns
True if an dataset of the given name with given dimension exists.

Definition at line 467 of file parallelio.h.

467 {
468 TRACE();
469
470 MBool returnValue = derived().b_hasDataset(name, dimension);
471
472 return returnValue;
473}

◆ hasObject()

template<class Backend >
MBool ParallelIoBase< Backend >::hasObject ( const MString path)

Definition at line 505 of file parallelio.h.

505 {
506 TRACE();
507
508 MBool returnValue = derived().b_hasObject(name);
509
510 return returnValue;
511}

◆ isParallelIoFile()

template<class Backend >
MBool ParallelIoBase< Backend >::isParallelIoFile ( const MString name,
const MPI_Comm &  mpiComm 
)
static
Author
Ramandeep Jain (HiWi) raman.nosp@m.deep.nosp@m.jain@.nosp@m.gmai.nosp@m.l.com, Konstantin Froehlich
Date
2015-03-01
Parameters
[in]mpiCommMPI communicator to be used.

This must be a collective call from all ranks in mpiComm.

Definition at line 364 of file parallelio.h.

364 {
365 TRACE();
366
367 return Backend::b_isValidFile(name, mpiComm);
368}

◆ readArray() [1/5]

template<class Backend >
template<class T >
void ParallelIoBase< Backend >::readArray ( std::vector< T > &  array,
const MString name,
size_type  memoryStride = -1,
size_type  diskStride = -1 
)
inline

Definition at line 195 of file parallelio.h.

195 {
196 readArray(array.data(), name, memoryStride, diskStride);
197 }
void readArray(T *array, const MString &name, size_type memoryStride=-1, size_type diskStride=-1)
Read array data from file. [MPI]
Definition: parallelio.h:1523

◆ readArray() [2/5]

template<class Backend >
template<class T >
void ParallelIoBase< Backend >::readArray ( T *  array,
const MString datasetName,
const MString name,
const size_type  noDims,
const size_type localCount 
)
Author
Marian Albers m.alb.nosp@m.ers@.nosp@m.aia.r.nosp@m.wth-.nosp@m.aache.nosp@m.n.de
Date
2022-09-029
Template Parameters
TData type of the array (can be either MFloat or MInt and must match the previous definition of the array).
Parameters
[in]arrayPointer to the data in memory.
[in]nameName of the array in the file.
[in]pathGroup name in which the array shall be located (works only with HDF5).
[in]noDimsNumber of dimensions of the array
[in]localCountSize of the array in each dimension, must of of size noDims

Definition at line 1645 of file parallelio.h.

1646 {
1647 TRACE();
1648
1649#ifdef MAIA_EXTRA_DEBUG
1650 // In debug mode, check if array was already read and give a warning
1651 if(m_readArrays.count(name) != 0) {
1652 std::cerr << "Warning: in " << AT_ << " (" << __FILE__ << ":" << __LINE__ << ") "
1653 << "the array '" << name << "' is read more than once. Make sure that this is the"
1654 << "intended behavior." << std::endl;
1655 }
1656 m_readArrays.insert(name);
1657#endif
1658
1659 if(array == 0 && localCount[0] > 0) {
1660 mTerm(1, AT_, "Data pointer must point to a valid location (is: null pointer)!");
1661 }
1662
1663 ScratchSpace<size_type> offset(noDims, FUN_, "offset");
1664 offset.fill(0);
1665
1666 // Read array
1667 derived().b_readArray(array, datasetName, name, noDims, &offset[0], localCount);
1668}
std::set< MString > m_readArrays
Definition: parallelio.h:258
This class is a ScratchSpace.
Definition: scratch.h:758

◆ readArray() [3/5]

template<class Backend >
template<class T >
void ParallelIoBase< Backend >::readArray ( T *  array,
const MString datasetName,
const MString name,
const size_type  noDims,
const size_type offset,
const size_type localCount 
)
Author
Marian Albers m.alb.nosp@m.ers@.nosp@m.aia.r.nosp@m.wth-.nosp@m.aache.nosp@m.n.de
Date
2022-09-029
Template Parameters
TData type of the array (can be either MFloat or MInt and must match the previous definition of the array).
Parameters
[in]arrayPointer to the data in memory.
[in]nameName of the array in the file.
[in]pathGroup name in which the array shall be located (works only with HDF5).
[in]noDimsNumber of dimensions of the array
[in]offsetOffset of the of the array in the dataset in each dimension, must be of size noDims
[in]localCountSize of the array in each dimension, must of of size noDims

Definition at line 1602 of file parallelio.h.

1607 {
1608 TRACE();
1609
1610#ifdef MAIA_EXTRA_DEBUG
1611 // In debug mode, check if array was already read and give a warning
1612 if(m_readArrays.count(name) != 0) {
1613 std::cerr << "Warning: in " << AT_ << " (" << __FILE__ << ":" << __LINE__ << ") "
1614 << "the array '" << name << "' is read more than once. Make sure that this is the"
1615 << "intended behavior." << std::endl;
1616 }
1617 m_readArrays.insert(name);
1618#endif
1619
1620 if(array == 0 && localCount[0] > 0) {
1621 mTerm(1, AT_, "Data pointer must point to a valid location (is: null pointer)!");
1622 }
1623
1624 // Read array
1625 derived().b_readArray(array, datasetName, name, noDims, offset, localCount);
1626}

◆ readArray() [4/5]

template<class Backend >
template<class T >
void ParallelIoBase< Backend >::readArray ( T *  array,
const MString name,
size_type  memoryStride = -1,
size_type  diskStride = -1 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-09
Template Parameters
TData type of the array (can be either MFloat or MInt and must match the previous definition of the array).
Parameters
[in]arrayPointer to the data in memory.
[in]nameName of the array in the file.
[in]memoryStrideOffset between two data values in memory (or -1 to set automatically).
[in]diskStrideOffset 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.

1523 {
1524 TRACE();
1525
1526 // Abort if offsets were not set
1527 if(!m_offsetsAreSet) {
1528 mTerm(1, AT_, "Offsets have to be set before calling readArray!");
1529 }
1530
1531#ifdef MAIA_EXTRA_DEBUG
1532 // In debug mode, check if array was already read and give a warning
1533 if(m_readArrays.count(name) != 0) {
1534 std::cerr << "Warning: in " << AT_ << " (" << __FILE__ << ":" << __LINE__ << ") "
1535 << "the array '" << name << "' is read more than once. Make sure that this is the"
1536 << "intended behavior." << std::endl;
1537 }
1538 m_readArrays.insert(name);
1539#endif
1540
1541 // Ensure reasonable parameter values
1542 if(memoryStride != -1 && memoryStride <= 0) {
1543 mTerm(1, AT_, "memoryStride must be greater than zero (or -1 to set automatically)!");
1544 }
1545 if(diskStride != -1 && diskStride <= 0) {
1546 mTerm(1, AT_, "diskStride must be greater than zero (or -1 to set automatically)!");
1547 }
1548 if(array == 0 && m_localCount[0] > 0) {
1549 mTerm(1, AT_, "Data pointer must point to a valid location (is: null pointer)!");
1550 }
1551
1552 // Set actual stride if no explicit stride value was given
1553 size_type actualMemoryStride;
1554 if(memoryStride == -1) {
1555 actualMemoryStride = 1;
1556 } else {
1557 actualMemoryStride = memoryStride;
1558 }
1559 size_type actualDiskStride = (diskStride == -1) ? 1 : diskStride;
1560
1561 // Set number of chunks
1562 const size_type noChunks = (m_noChunks == -1) ? 1 : m_noChunks;
1563
1564 size_type noDims = getDatasetNoDims(name);
1565
1566 ScratchSpace<MPI_Offset> localCount(noDims, FUN_, "localCount");
1567 ScratchSpace<MPI_Offset> offset(noDims, FUN_, "offset");
1568
1569 localCount[0] = m_localCount[0];
1570 // Set counts for all dimensions except first to total dimension length
1571 // m_localCount[0] is set via setOffset(...)
1572 for(MInt dimId = 1; dimId < noDims; dimId++) {
1573 localCount[dimId] = getArraySize(name, dimId);
1574 }
1575 for(MInt dimId = 0; dimId < noDims; dimId++) {
1576 offset[dimId] = m_offset[dimId];
1577 }
1578
1579 // Read array
1580 derived().b_readArray(array, name, noDims, &offset[0], &localCount[0], actualMemoryStride, noChunks,
1581 actualDiskStride);
1582}
size_type m_noChunks
Definition: parallelio.h:247
MBool m_offsetsAreSet
Definition: parallelio.h:244
std::vector< size_type > m_offset
Definition: parallelio.h:246
std::vector< size_type > m_localCount
Definition: parallelio.h:245

◆ readArray() [5/5]

template<class Backend >
template<class T >
void ParallelIoBase< Backend >::readArray ( T *  array,
const std::vector< MString > &  names,
size_type  memoryStride = -1,
size_type  diskStride = -1 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-09
Template Parameters
TData type of the arrays (can be either MFloat or MInt and must match the previous definition of the array).
Parameters
[in]arrayPointer to the data of the first array in memory.
[in]namesNames of the arrays in the file.
[in]memoryStrideOffset between two data values in memory (or -1 to set automatically).
[in]diskStrideOffset 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.

1712 {
1713 TRACE();
1714
1715 // Ensure reasonable parameter values
1716 if(memoryStride != -1 && memoryStride < static_cast<size_type>(names.size())) {
1717 mTerm(1, AT_,
1718 "memoryStride must be greater than or equal to number of "
1719 "variables (or -1 to set automatically)!");
1720 }
1721 if(diskStride != -1 && diskStride < static_cast<size_type>(names.size())) {
1722 mTerm(1, AT_,
1723 "diskStride must be greater than or equal to number of "
1724 "variables (or -1 to set automatically)!");
1725 }
1726
1727 // Set actual stride if no explicit stride value was given
1728 size_type actualMemoryStride;
1729 if(memoryStride == -1) {
1730 actualMemoryStride = static_cast<size_type>(names.size());
1731 } else {
1732 actualMemoryStride = memoryStride;
1733 }
1734 size_type actualDiskStride = (diskStride == -1) ? 1 : diskStride;
1735
1736 // Read arrays from vector using normal readArray call
1737 for(std::vector<MString>::size_type i = 0; i < names.size(); i++) {
1738 readArray(array + i, names[i], actualMemoryStride, actualDiskStride);
1739 }
1740}

◆ readScalar()

template<class Backend >
template<class T >
void ParallelIoBase< Backend >::readScalar ( T *  scalar,
const MString name 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-09
Template Parameters
TData type of the scalar (can be either MFloat or MInt and must match the previous definition of the scalar).
Parameters
[in]scalarValue that should be read from file.
[in]nameName of the scalar in the file.

Note: This method is the exact counterpart to writeScalar().

Definition at line 1758 of file parallelio.h.

1758 {
1759#ifdef MAIA_EXTRA_DEBUG
1760 // In debug mode, check if scalar was already read and give a warning
1761 if(m_readScalars.count(name) != 0) {
1762 std::cerr << "Warning: in " << AT_ << " (" << __FILE__ << ":" << __LINE__ << ") "
1763 << "the scalar '" << name << "' is read more than once. Make sure that this is the"
1764 << "intended behavior." << std::endl;
1765 }
1766 m_readScalars.insert(name);
1767#endif
1768
1769 derived().b_readScalar(scalar, name);
1770}
std::set< MString > m_readScalars
Definition: parallelio.h:259

◆ setAttribute() [1/2]

template<class Backend >
void ParallelIoBase< Backend >::setAttribute ( const MChar value,
const MString name,
const MString datasetName = "" 
)

Definition at line 1445 of file parallelio.h.

1445 {
1446 const MString tmpValue = MString(value);
1447 setAttributes(&tmpValue, name, 1, datasetName);
1448}
void setAttributes(const T *value, const MString &name, size_type totalCount, const MString &datasetName="")
Definition: parallelio.h:1452
std::basic_string< char > MString
Definition: maiatypes.h:55

◆ setAttribute() [2/2]

template<class Backend >
template<class T >
void ParallelIoBase< Backend >::setAttribute ( const T &  value,
const MString name,
const MString datasetName = "" 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-09
Template Parameters
TData type of the attribute (can be either MFloat, MInt, or MString).
Parameters
[in]valueValue of the attribute.
[in]nameName of the attribute (must not be empty).
[in]datasetNameAttaches 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.

1438 {
1439 TRACE();
1440 setAttributes(&value, name, 1, datasetName);
1441}

◆ setAttributes()

template<class Backend >
template<class T >
void ParallelIoBase< Backend >::setAttributes ( const T *  value,
const MString name,
size_type  totalCount,
const MString datasetName = "" 
)

Definition at line 1452 of file parallelio.h.

1455 {
1456 TRACE();
1457
1458#ifdef DISABLE_OUTPUT
1459 return;
1460#endif
1461
1462 using namespace maia::parallel_io;
1463
1464 // Prevent rewriting data in read-only mode
1465 if(m_fileMode == PIO_READ) {
1466 mTerm(1, AT_, "Cannot write data in read-only mode!");
1467 }
1468
1469 // Ensure reasonable parameter values
1470 if(name.empty()) {
1471 mTerm(1, AT_, "Invalid name! Name must not be empty.");
1472 }
1473 ASSERT(totalCount > 0, "");
1474
1475 derived().b_setAttribute(value, name, datasetName, totalCount);
1476}

◆ setOffset() [1/3]

template<class Backend >
void ParallelIoBase< Backend >::setOffset ( const size_type  localCount,
const size_type  offset,
const size_type  noDims 
)

Definition at line 1889 of file parallelio.h.

1889 {
1890 TRACE();
1891
1892 setOffset(localCount, offset, noDims, -1);
1893}
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.
Definition: parallelio.h:1840

◆ setOffset() [2/3]

template<class Backend >
void ParallelIoBase< Backend >::setOffset ( const size_type  localCount,
const size_type  offset,
const size_type  noDims,
const size_type  noChunks 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-09
Parameters
[in]localCountNumber of elements in an array on the current domain.
[in]offsetStarting offset in the file, from which this domain should start writing/reading.
[in]noDimsnumber 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]noChunksNumber 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.

1843 {
1844 TRACE();
1845
1846 // Ensure reasonable parameter values
1847 if(localCount < 0) {
1848 TERMM(1, "The local size may not be less than zero (was:" + std::to_string(localCount) + ")!");
1849 }
1850 if(offset < 0) {
1851 TERMM(1, "The offset may not be less than zero (was:" + std::to_string(offset) + ")!");
1852 }
1853 if(noDims <= 0) {
1854 TERMM(1, "The number of dimensions must at least be one!");
1855 }
1856 if(noChunks == 0) {
1857 TERMM(1, "The number of chunks may not be zero!");
1858 }
1859
1860 // TODO labels:IO add support for setting a multi-dimensional count/offset for reading/writing data
1861
1862 // Set member variables
1863 m_localCount.reserve(noDims);
1864 m_localCount.assign(noDims, -1); // -1 as a precaution
1865 m_localCount[0] = localCount;
1866
1867 // Note:
1868 // if noDims>1: {m_localCount[i], i>0} will be set to the total count for that
1869 // dimension before reading/writing the multi-d array
1870 // setting the local count for multiple dimensions is not yet supported
1871 m_offset.reserve(noDims);
1872 m_offset.assign(noDims, 0);
1873 m_offset[0] = offset;
1874
1875 if(noChunks > 0) {
1876 m_noChunks = noChunks;
1877 }
1878
1879 // Set state variable
1880 m_offsetsAreSet = true;
1881}

◆ setOffset() [3/3]

template<class Backend >
void ParallelIoBase< Backend >::setOffset ( size_type  localCount,
size_type  offset 
)

Definition at line 1901 of file parallelio.h.

1901 {
1902 TRACE();
1903
1904 setOffset(localCount, offset, 1, -1);
1905}

◆ validateType()

template<class Backend >
void ParallelIoBase< Backend >::validateType ( const maiabd_type  type) const
protected
Author
Michael Schlottke-Lakemper (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2017-06-23
Parameters
[in]typeThe type constant (integer) to check for validity.

Definition at line 1976 of file parallelio.h.

1976 {
1977 TRACE();
1978
1979 using namespace maia::parallel_io;
1980
1981 if(type != PIO_FLOAT && type != PIO_INT && type != PIO_LONG && type != PIO_STRING && type != PIO_UCHAR
1982 && type != PIO_ULONGLONG) {
1983 TERMM(1, "Invalid data type! Must be PIO_FLOAT or PIO_INT or PIO_LONG or PIO_STRING or PIO_UCHAR.");
1984 }
1985}

◆ writeArray() [1/4]

template<class Backend >
template<class T >
void ParallelIoBase< Backend >::writeArray ( const T *  array,
const MString name,
size_type  memoryStride = -1,
size_type  diskStride = -1 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-09
Template Parameters
TData type of the array (can be either MFloat or MInt and must match the previous definition of the array).
Parameters
[in]arrayPointer to the data in memory.
[in]nameName of the array in the file.
[in]memoryStrideOffset between two data values in memory (or -1 to set automatically).
[in]diskStrideOffset 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.

1047 {
1048 TRACE();
1049
1050#ifdef DISABLE_OUTPUT
1051 return;
1052#endif
1053
1054 using namespace maia::parallel_io;
1055
1056 // Prevent rewriting data in read-only mode
1057 if(m_fileMode == PIO_READ) {
1058 TERMM(1, "Cannot write data in read-only mode!");
1059 }
1060
1061 // Abort if offsets were not set
1062 if(!m_offsetsAreSet) {
1063 TERMM(1, "Offsets have to be set before calling writeArray!");
1064 }
1065
1066#ifdef MAIA_EXTRA_DEBUG
1067 // In debug mode, check if array was already written and give a warning
1068 if(m_writtenArrays.count(name) != 0) {
1069 std::cerr << "Warning: in " << AT_ << " (" << __FILE__ << ":" << __LINE__ << ") "
1070 << "the array '" << name << "' is written more than once. Make sure that this is the"
1071 << "intended behavior." << std::endl;
1072 }
1073 m_writtenArrays.insert(name);
1074 if(m_unwrittenArrays.count(name)) {
1075 m_unwrittenArrays.erase(name);
1076 }
1077#endif
1078
1079 // Ensure reasonable parameter values
1080 if(memoryStride != -1 && memoryStride <= 0) {
1081 TERMM(1, "memoryStride must be greater than zero (or -1 to set automatically)!");
1082 }
1083 if(diskStride != -1 && diskStride <= 0) {
1084 TERMM(1, "diskStride must be greater than zero (or -1 to set automatically)!");
1085 }
1086 if(array == 0 && m_localCount[0] > 0) {
1087 TERMM(1, "Data pointer must point to a valid location (is: null pointer)!");
1088 }
1089
1090 // Set actual stride if no explicit stride value was given
1091 size_type actualMemoryStride;
1092 if(memoryStride == -1) {
1093 actualMemoryStride = 1;
1094 } else {
1095 actualMemoryStride = memoryStride;
1096 }
1097 size_type actualDiskStride = (diskStride == -1) ? 1 : diskStride;
1098
1099 // Save header (if not yet
1100 if(!m_isHeaderSaved) {
1101 derived().b_saveHeader();
1102 m_isHeaderSaved = true;
1103 }
1104
1105 // Set number of chunks
1106 const size_type noChunks = (m_noChunks == -1) ? 1 : m_noChunks;
1107
1108 size_type noDims = getDatasetNoDims(name);
1109
1110 ScratchSpace<MPI_Offset> localCount(noDims, FUN_, "localCount");
1111 ScratchSpace<MPI_Offset> offset(noDims, FUN_, "offset");
1112
1113 // Set counts for all dimensions except first to total dimension length
1114 // localCount[0] is set via setOffset(...)
1115 localCount[0] = m_localCount[0];
1116 for(MInt dimId = 1; dimId < noDims; dimId++) {
1117 localCount[dimId] = getArraySize(name, dimId);
1118 }
1119 for(MInt dimId = 0; dimId < noDims; dimId++) {
1120 offset[dimId] = m_offset[dimId];
1121 }
1122
1123 // Write array
1124 derived().b_writeArray(array, name, noDims, &offset[0], &localCount[0], actualMemoryStride, noChunks,
1125 actualDiskStride);
1126}
std::set< MString > m_writtenArrays
Definition: parallelio.h:256

◆ writeArray() [2/4]

template<class Backend >
template<class T >
void ParallelIoBase< Backend >::writeArray ( const T *  array,
const MString path,
const MString name,
const size_type  noDims,
size_type offset,
size_type localCount 
)
Author
Marian Albers m.alb.nosp@m.ers@.nosp@m.aia.r.nosp@m.wth-.nosp@m.aache.nosp@m.n.de
Date
2022-09-029
Template Parameters
TData type of the array (can be either MFloat or MInt and must match the previous definition of the array).
Parameters
[in]arrayPointer to the data in memory.
[in]pathGroup name in which the array shall be located (works only with HDF5).
[in]nameName of the array in the file.
[in]noDimsNumber of dimensions of the array
[in]offsetOffset of the of the array in the dataset in each dimension, must be of size noDims
[in]localCountSize of the array in each dimension, must of of size noDims

Definition at line 1255 of file parallelio.h.

1256 {
1257 TRACE();
1258
1259 ScratchSpace<size_type> ghost(noDims, FUN_, "ghost");
1260 ghost.fill(0);
1261
1262 derived().b_writeArray(array, path, name, noDims, offset, localCount, &ghost[0]);
1263}

◆ writeArray() [3/4]

template<class Backend >
template<class T >
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 
)
Author
Marian Albers m.alb.nosp@m.ers@.nosp@m.aia.r.nosp@m.wth-.nosp@m.aache.nosp@m.n.de
Date
2022-09-029
Template Parameters
TData type of the array (can be either MFloat or MInt and must match the previous definition of the array).
Parameters
[in]arrayPointer to the data in memory.
[in]pathGroup name in which the array shall be located (works only with HDF5).
[in]nameName of the array in the file.
[in]noDimsNumber of dimensions of the array
[in]offsetOffset of the of the array in the dataset in each dimension, must be of size noDims
[in]localCountSize of the array in each dimension, must of of size noDims
[in]ghostNumber 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.

1298 {
1299 TRACE();
1300
1301#ifdef DISABLE_OUTPUT
1302 return;
1303#endif
1304
1305 using namespace maia::parallel_io;
1306
1307 // Prevent rewriting data in read-only mode
1308 if(m_fileMode == PIO_READ) {
1309 TERMM(1, "Cannot write data in read-only mode!");
1310 }
1311
1312#ifdef MAIA_EXTRA_DEBUG
1313 // In debug mode, check if array was already written and give a warning
1314 if(m_writtenArrays.count(name) != 0) {
1315 std::cerr << "Warning: in " << AT_ << " (" << __FILE__ << ":" << __LINE__ << ") "
1316 << "the array '" << name << "' is written more than once. Make sure that this is the"
1317 << "intended behavior." << std::endl;
1318 }
1319 m_writtenArrays.insert(name);
1320 if(m_unwrittenArrays.count(name)) {
1321 m_unwrittenArrays.erase(name);
1322 }
1323#endif
1324
1325 // Save header (if not yet
1326 if(!m_isHeaderSaved) {
1327 derived().b_saveHeader();
1328 m_isHeaderSaved = true;
1329 }
1330
1331 // Ensure reasonable parameter values
1332 if(array == 0 && localCount[0] > 0) {
1333 TERMM(1, "Data pointer must point to a valid location (is: null pointer)!");
1334 }
1335
1336 // Write array
1337 derived().b_writeArray(array, path, name, noDims, offset, localCount, ghost);
1338}

◆ writeArray() [4/4]

template<class Backend >
template<class T >
void ParallelIoBase< Backend >::writeArray ( const T *  array,
const std::vector< MString > &  names,
size_type  memoryStride = -1,
size_type  diskStride = -1 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-09
Template Parameters
TData type of the arrays (can be either MFloat or MInt and must match the previous definition of the array).
Parameters
[in]arrayPointer to the data of the first array in memory.
[in]namesNames of the arrays in the file.
[in]memoryStrideOffset between two data values in memory (or -1 to set automatically).
[in]diskStrideOffset 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.

1205 {
1206 TRACE();
1207
1208#ifdef DISABLE_OUTPUT
1209 return;
1210#endif
1211
1212 // Ensure reasonable parameter values
1213 if(memoryStride != -1 && memoryStride < static_cast<size_type>(names.size())) {
1214 TERMM(1, "memoryStride must be greater than or equal to number of "
1215 "variables (or -1 to set automatically)!");
1216 }
1217 if(diskStride != -1 && diskStride < static_cast<size_type>(names.size())) {
1218 TERMM(1, "diskStride must be greater than or equal to number of "
1219 "variables (or -1 to set automatically)!");
1220 }
1221
1222 // Set actual stride if no explicit stride value was given
1223 size_type actualMemoryStride;
1224 if(memoryStride == -1) {
1225 actualMemoryStride = static_cast<size_type>(names.size());
1226 } else {
1227 actualMemoryStride = memoryStride;
1228 }
1229 size_type actualDiskStride = (diskStride == -1) ? 1 : diskStride;
1230
1231 // Write arrays from vector using normal writeArray call
1232 for(std::vector<MString>::size_type i = 0; i < names.size(); i++) {
1233 writeArray(array + i, names[i], actualMemoryStride, actualDiskStride);
1234 }
1235}
void writeArray(const T *array, const MString &name, size_type memoryStride=-1, size_type diskStride=-1)
Write array data to file. [MPI]
Definition: parallelio.h:1046

◆ writeScalar()

template<class Backend >
template<class T >
void ParallelIoBase< Backend >::writeScalar ( scalar,
const MString name 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de, Konstantin Froehlich
Date
2013-04-09
Template Parameters
TData type of the scalar (can be either MFloat or MInt and must match the previous definition of the scalar).
Parameters
[in]scalarValue that should be written to file.
[in]nameName of the scalar in the file.

Definition at line 1354 of file parallelio.h.

1354 {
1355 TRACE();
1356
1357#ifdef DISABLE_OUTPUT
1358 return;
1359#endif
1360
1361 using namespace maia::parallel_io;
1362
1363 // Prevent rewriting data in read-only mode
1364 if(m_fileMode == PIO_READ) {
1365 mTerm(1, AT_, "Cannot write data in read-only mode!");
1366 }
1367
1368#ifndef NDEBUG
1369 // In debug mode check if scalar value is the same on all ranks
1370 const T localValue = scalar;
1371 ScratchSpace<T> globalValues(m_noDomains, AT_, "globalValues");
1372 // Gather scalar values on rank 0
1373 MPI_Gather(&localValue, 1, maia::type_traits<T>::mpiType(), &globalValues[0], 1, maia::type_traits<T>::mpiType(), 0,
1374 m_mpiComm, AT_, "localValue", "globalValues[0]");
1375
1376 // Compare values
1377 if(m_domainId == 0) {
1378 for(MInt i = 0; i < m_noDomains; i++) {
1379 MBool equal = true;
1380 // Note: cannot use == for floating point values
1381 if(globalValues[i] < localValue) {
1382 equal = false;
1383 }
1384 if(globalValues[i] > localValue) {
1385 equal = false;
1386 }
1387 if(!equal) {
1388 TERMM(1, "Error: scalar variable has not the same value on all domains, rank 0: " + std::to_string(localValue)
1389 + " != " + std::to_string(globalValues[i]) + " (domain " + std::to_string(i) + ")");
1390 }
1391 }
1392 }
1393#endif
1394
1395#ifdef MAIA_EXTRA_DEBUG
1396 // In extra debug mode, check if scalar was already written and give a warning
1397 if(m_writtenScalars.count(name) != 0) {
1398 std::cerr << "Warning: in " << AT_ << " (" << __FILE__ << ":" << __LINE__ << ") "
1399 << "the scalar '" << name << "' is written more than once. Make sure that this is the"
1400 << "intended behavior." << std::endl;
1401 }
1402 m_writtenScalars.insert(name);
1403 if(m_unwrittenScalars.count(name)) {
1404 m_unwrittenScalars.erase(name);
1405 }
1406#endif
1407
1408 // Save header (if not yet done)
1409 if(!m_isHeaderSaved) {
1410 derived().b_saveHeader();
1411 m_isHeaderSaved = true;
1412 }
1413
1414 derived().b_writeScalar(scalar, name);
1415}
std::set< MString > m_writtenScalars
Definition: parallelio.h:257
int MPI_Gather(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, int recvcount, MPI_Datatype recvtype, int root, MPI_Comm comm, const MString &name, const MString &sndvarname, const MString &rcvvarname)
same as MPI_Gather

Member Data Documentation

◆ m_creationTime

template<class Backend >
MFloat ParallelIoBase< Backend >::m_creationTime = MFloatNaN
protected

Definition at line 250 of file parallelio.h.

◆ m_domainId

template<class Backend >
MInt ParallelIoBase< Backend >::m_domainId = -1
protected

Definition at line 236 of file parallelio.h.

◆ m_fileMode

template<class Backend >
const MInt ParallelIoBase< Backend >::m_fileMode = -1
protected

Definition at line 231 of file parallelio.h.

◆ m_fileName

template<class Backend >
const MString ParallelIoBase< Backend >::m_fileName {}
protected

Definition at line 230 of file parallelio.h.

◆ m_isHeaderSaved

template<class Backend >
MBool ParallelIoBase< Backend >::m_isHeaderSaved = false
protected

Definition at line 232 of file parallelio.h.

◆ m_localCount

template<class Backend >
std::vector<size_type> ParallelIoBase< Backend >::m_localCount = {-1}
protected

Definition at line 245 of file parallelio.h.

◆ m_mpiComm

template<class Backend >
MPI_Comm ParallelIoBase< Backend >::m_mpiComm {}
protected

Definition at line 240 of file parallelio.h.

◆ m_mpiInfo

template<class Backend >
MPI_Info ParallelIoBase< Backend >::m_mpiInfo = MPI_INFO_NULL
protected

Definition at line 241 of file parallelio.h.

◆ m_noChunks

template<class Backend >
size_type ParallelIoBase< Backend >::m_noChunks = -1
protected

Definition at line 247 of file parallelio.h.

◆ m_noDomains

template<class Backend >
MInt ParallelIoBase< Backend >::m_noDomains = -1
protected

Definition at line 237 of file parallelio.h.

◆ m_offset

template<class Backend >
std::vector<size_type> ParallelIoBase< Backend >::m_offset = {-1}
protected

Definition at line 246 of file parallelio.h.

◆ m_offsetsAreSet

template<class Backend >
MBool ParallelIoBase< Backend >::m_offsetsAreSet = false
protected

Definition at line 244 of file parallelio.h.

◆ m_readArrays

template<class Backend >
std::set<MString> ParallelIoBase< Backend >::m_readArrays
protected

Definition at line 258 of file parallelio.h.

◆ m_readScalars

template<class Backend >
std::set<MString> ParallelIoBase< Backend >::m_readScalars
protected

Definition at line 259 of file parallelio.h.

◆ m_unwrittenArrays

template<class Backend >
std::set<MString> ParallelIoBase< Backend >::m_unwrittenArrays
protected

Definition at line 254 of file parallelio.h.

◆ m_unwrittenScalars

template<class Backend >
std::set<MString> ParallelIoBase< Backend >::m_unwrittenScalars
protected

Definition at line 255 of file parallelio.h.

◆ m_writtenArrays

template<class Backend >
std::set<MString> ParallelIoBase< Backend >::m_writtenArrays
protected

Definition at line 256 of file parallelio.h.

◆ m_writtenScalars

template<class Backend >
std::set<MString> ParallelIoBase< Backend >::m_writtenScalars
protected

Definition at line 257 of file parallelio.h.

◆ s_maxChunkSize

template<class Backend >
const size_type ParallelIoBase< Backend >::s_maxChunkSize = 2147483645
staticprotected

Definition at line 265 of file parallelio.h.

◆ s_mpiSizeType

template<typename Backend >
const MPI_Datatype ParallelIoBase< Backend >::s_mpiSizeType = maia::type_traits<ParallelIoBase<Backend>::size_type>::mpiType()
staticprotected

Definition at line 270 of file parallelio.h.


The documentation for this class was generated from the following file: