MAIA bb96820c
Multiphysics at AIA
Loading...
Searching...
No Matches
ParallelIoPNetcdf Class Reference

#include <parallelio_pnetcdf.h>

Inheritance diagram for ParallelIoPNetcdf:
[legend]
Collaboration diagram for ParallelIoPNetcdf:
[legend]

Classes

struct  NcDimProxy
 

Public Member Functions

 ParallelIoPNetcdf (const MString &fileName, MInt fileMode, const MPI_Comm &mpiComm)
 Creates a new object to read and write big data files. [MPI] More...
 
 ~ParallelIoPNetcdf () override
 Calls close(). [MPI] More...
 
void close ()
 Close the file (normally called by the destructor but needs to be explicitely called earlier in specific cases, e.g. when. More...
 
template<>
void b_setAttribute (const MString *value, const MString &name, const MString &datasetName, const size_type totalCount)
 Set an attribute in the file (string version). More...
 
template<>
void b_writeArray (const MString *array, const MString &name, const size_type noDims, MPI_Offset *start, MPI_Offset *count, MPI_Offset memoryStride, const size_type noChunks, MPI_Offset diskStride)
 Writes array data to file (string version). [MPI] More...
 
template<>
void b_writeScalar (const MString &scalar, const MString &name)
 Writes scalar data to file (string version). [MPI] More...
 
template<>
void b_readArray (MString *array, const MString &name, const size_type noDims, MPI_Offset *start, MPI_Offset *count, MPI_Offset memoryStride, const size_type noChunks, MPI_Offset diskStride)
 Read array data from file (string version). [MPI] More...
 
template<>
void b_readScalar (MString *scalar, const MString &name)
 Read scalar data from file (float version). [MPI] More...
 
template<>
void b_getAttribute (MString *const value, const MString &name, const MString &datasetName, const size_type totalCount)
 Retrieve an attribute from file (string version). More...
 
- Public Member Functions inherited from ParallelIoBase< ParallelIoPNetcdf >
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...
 
void writeArray (const T *array, const MString &name, size_type memoryStride=-1, size_type diskStride=-1)
 Write array data to file. [MPI] More...
 
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...
 
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...
 
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...
 
void writeScalar (T scalar, const MString &name)
 Write scalar data to file. [MPI] More...
 
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="")
 
void setAttributes (const T *value, const MString &name, size_type totalCount, const MString &datasetName="")
 
void readArray (T *array, const MString &name, size_type memoryStride=-1, size_type diskStride=-1)
 Read array data from file. [MPI] More...
 
void readArray (std::vector< T > &array, const MString &name, size_type memoryStride=-1, size_type diskStride=-1)
 
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...
 
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...
 
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...
 
void readScalar (T *scalar, const MString &name)
 Read scalar data from file. [MPI] More...
 
void getAttribute (T *value, const MString &name, const MString &datasetName="")
 Retrieve a file or dataset attribute. More...
 
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...
 

Private Types

using NcDimMap = std::map< size_type, NcDimProxy >
 

Private Member Functions

void b_ncEndDef ()
 Leave the define mode (NetCDF-specific). [MPI] More...
 
void b_ncRedef ()
 Enter define mode (NetCDF-specific). [MPI] More...
 
void b_saveHeader ()
 Adds all additional header information that are needed in the file. More...
 
void b_addAdditionalHeader ()
 Write additional headers to file (e.g. grid file name, creation date etc.). [MPI] More...
 
void b_writeAdditionalData ()
 Write additional data to file (NetCDF-specific). [MPI] More...
 
void b_defineArray (maiabd_type type, const MString &name, size_type noDims, size_type *totalCount)
 Defines an array in the file. [MPI] More...
 
void b_defineArray (maiabd_type, const MString &, const MString &, const size_type, const size_type *)
 
void b_defineScalar (maiabd_type type, const MString &name)
 Defines a scalar in the file. [MPI] More...
 
MBool b_hasDataset (const MString &name, const size_type dimension)
 Check if the file contains a dataset with the given name and number of dimensions. More...
 
MBool b_hasDataset (const MString &name, const MString &path)
 
MBool b_hasObject (const MString &path)
 
MInt b_getDatasetType (const MString &name)
 Returns the data type of a dataset in the file (can be array, multi-D array or scalar). More...
 
void b_getDatasetNames (std::vector< MString > &names, const size_type dimension)
 Returns a vector of all existing datasets with the given number of dimensions in the file (if any). [MPI] More...
 
void b_getDatasetNames (std::vector< MString > &names, const MString &path)
 
void b_getGroupNames (std::vector< MString > &names, const MString &path)
 
size_type b_getDatasetNoDims (const MString &name)
 Get number of dimensions of a dataset with the given name. More...
 
size_type b_getDatasetNoDims (const MString &path, const MString &name)
 
size_type b_getDatasetSize (const MString &name, const size_type dimensionId=0)
 Get the length of one dimension of an arbitrary array in the file. More...
 
void b_getDatasetSize (const MString &name, const MString &path, size_type noDims, size_type *data)
 
MBool b_hasAttribute (const MString &name, const MString &datasetName="")
 Check if a given attribute exists in the file. More...
 
MInt b_getAttributeType (const MString &name, const MString &datasetName="")
 Returns the data type of an attribute in the file. More...
 
void b_printFileHints ()
 Print PNetCDF file hints to cerr. More...
 
template<class T >
void b_writeArray (const T *const array, const MString &name, const size_type noDims, MPI_Offset *start, MPI_Offset *count, MPI_Offset memoryStride, const size_type noChunks, MPI_Offset diskStride)
 Writes array data to file (generic version). [MPI] More...
 
template<class T >
void b_writeArray (const T *, const MString &, const MString &, const size_type, const size_type *, const size_type *, const size_type *)
 
template<class T >
void b_writeScalar (const T scalar, const MString &name)
 Writes scalar data to file (generic version). [MPI] More...
 
template<class T >
void b_readArray (T *array, const MString &name, const size_type noDims, MPI_Offset *start, MPI_Offset *count, MPI_Offset memoryStride, const size_type noChunks, MPI_Offset diskStride)
 Read array data from file (generic version). [MPI] More...
 
template<class T >
void b_readArray (T *, const MString &, const MString &, const size_type, const size_type *, const size_type *)
 
template<class T >
void b_readScalar (T *scalar, const MString &name)
 Read scalar data from file (generic version). [MPI] More...
 
template<class T >
void b_setAttribute (const T *value, const MString &name, const MString &datasetName="", const size_type totalCount=1)
 Set an attribute in the file (generic version). More...
 
template<class T >
void b_getAttribute (T *const value, const MString &name, const MString &datasetName="", const size_type totalCount=1)
 Retrieve an attribute from file (generic version). More...
 
void b_getAttributeCount (const MString &name, size_type *totalCount, const MString &datasetName="")
 

Static Private Member Functions

static MBool b_isValidFile (const MString &name, const MPI_Comm &mpiComm)
 Check if specified file is a valid HDF5 file (i.e. can be opened). [MPI] More...
 
static const MString b_fileExt ()
 Returns backend-specific ending of filename (either ".Netcdf" or ".Hdf5") More...
 
static void b_error (MInt status, const MString &name, const MString &location)
 Check the status code of a HDF5 operation and output a meaningful message. More...
 

Private Attributes

friend ParallelIoBase< ParallelIoPNetcdf >
 
MInt b_ncId = -1
 
MBool b_ncDataMode = false
 
NcDimMap b_ncDimensions {}
 

Additional Inherited Members

- Public Types inherited from ParallelIoBase< ParallelIoPNetcdf >
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...
 
- Static Public Member Functions inherited from ParallelIoBase< ParallelIoPNetcdf >
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 inherited from ParallelIoBase< ParallelIoPNetcdf >
 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 inherited from ParallelIoBase< ParallelIoPNetcdf >
const MString m_fileName
 
const MInt m_fileMode
 
MBool m_isHeaderSaved
 
MInt m_domainId
 
MInt m_noDomains
 
MPI_Comm m_mpiComm
 
MPI_Info m_mpiInfo
 
MBool m_offsetsAreSet
 
std::vector< size_typem_localCount
 
std::vector< size_typem_offset
 
size_type m_noChunks
 
MFloat m_creationTime
 
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 inherited from ParallelIoBase< ParallelIoPNetcdf >
static const size_type s_maxChunkSize
 
static const MPI_Datatype s_mpiSizeType
 

Detailed Description

Definition at line 12 of file parallelio_pnetcdf.h.

Member Typedef Documentation

◆ NcDimMap

using ParallelIoPNetcdf::NcDimMap = std::map<size_type, NcDimProxy>
private

Definition at line 96 of file parallelio_pnetcdf.h.

Constructor & Destructor Documentation

◆ ParallelIoPNetcdf()

ParallelIoPNetcdf::ParallelIoPNetcdf ( const MString fileName,
MInt  fileMode,
const MPI_Comm &  mpiComm 
)
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2013-03-2
Author
(Modified by: )Ramandeep Jain (HiWi) raman.nosp@m.deep.nosp@m.jain@.nosp@m.gmai.nosp@m.l.com
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 372 of file parallelio_pnetcdf.cpp.

373 : ParallelIoBase<ParallelIoPNetcdf>(fileName, fileMode, mpiComm) {
374 TRACE();
375
376#ifdef DISABLE_OUTPUT
377 if(m_fileMode != PIO_READ) return;
378#endif
379
380 switch(m_fileMode) {
381 case PIO_CREATE: {
382 // Create a new file (do not overwrite existing)
383 MInt status =
384 ncmpi_create(m_mpiComm, m_fileName.c_str(), NC_NOCLOBBER | MAIA_NCMPI_FILE_TYPE, globalMpiInfo(), &b_ncId);
385 b_error(status, m_fileName, AT_);
386 } break;
387
388 case PIO_REPLACE: {
389 // Create a new file (overwrite existing)
390 MInt status =
391 ncmpi_create(m_mpiComm, m_fileName.c_str(), NC_CLOBBER | MAIA_NCMPI_FILE_TYPE, globalMpiInfo(), &b_ncId);
392 b_error(status, m_fileName, AT_);
393 } break;
394
395 case PIO_APPEND: {
396 // Attempt to open an existing file to append data
397 MInt status =
398 ncmpi_open(m_mpiComm, m_fileName.c_str(), NC_WRITE | MAIA_NCMPI_FILE_TYPE, globalMpiInfo(), &b_ncId);
399 b_error(status, m_fileName, AT_);
400
401 // Set correct data mode status
402 b_ncDataMode = true;
403
404 // Read in all existing dimensions and populate dimensions map
405 int noDims;
406 status = ncmpi_inq_ndims(b_ncId, &noDims);
407 b_error(status, m_fileName, AT_);
408 for(int dimId = 0; dimId < noDims; dimId++) {
409 MPI_Offset dimLength;
410 status = ncmpi_inq_dimlen(b_ncId, dimId, &dimLength);
411 b_error(status, m_fileName, AT_);
412 if(b_ncDimensions.count(dimLength) == 0u) {
413 b_ncDimensions[dimLength] = NcDimProxy();
414 b_ncDimensions[dimLength].dimNo = dimLength;
415 b_ncDimensions[dimLength].dimId = dimId;
416 }
417 }
418
419 // Enter definition mode so that new attributes/variables may be added
420 b_ncRedef();
421 } break;
422
423 case PIO_READ: {
424 // Open file for reading
425 MInt status = ncmpi_open(m_mpiComm, m_fileName.c_str(), NC_NOWRITE, globalMpiInfo(), &b_ncId);
426 b_error(status, m_fileName, AT_);
427
428 // Set correct data mode status
429 b_ncDataMode = true;
430 } break;
431
432 default: {
433 mTerm(1, AT_, "Unsupported file mode.");
434 } break;
435 }
436
437#ifdef MAIA_NCMPI_PRINT_FILE_HINTS
438 if(m_domainId == 0) {
439 std::cerr << std::endl << "Created/replaced/opened file: " << m_fileName << std::endl;
441 }
442#endif
443}
void b_printFileHints()
Print PNetCDF file hints to cerr.
friend ParallelIoBase< ParallelIoPNetcdf >
static void b_error(MInt status, const MString &name, const MString &location)
Check the status code of a HDF5 operation and output a meaningful message.
void b_ncRedef()
Enter define mode (NetCDF-specific). [MPI]
void mTerm(const MInt errorCode, const MString &location, const MString &message)
Definition: functions.cpp:29
const MPI_Info & globalMpiInfo()
Return global MPI information.
int32_t MInt
Definition: maiatypes.h:62
const MInt PIO_REPLACE
Definition: parallelio.h:36
const MInt PIO_APPEND
Definition: parallelio.h:38
const MInt PIO_CREATE
< File mode to create a new file. Aborts if file already exists
Definition: parallelio.h:34
const MInt PIO_READ
Definition: parallelio.h:40

◆ ~ParallelIoPNetcdf()

ParallelIoPNetcdf::~ParallelIoPNetcdf ( )
override
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

Definition at line 452 of file parallelio_pnetcdf.cpp.

452 {
453 TRACE();
454
455 // Check if file is already closed
456 if(b_ncId != -1) {
457 close();
458 }
459}
void close()
Close the file (normally called by the destructor but needs to be explicitely called earlier in speci...

Member Function Documentation

◆ b_addAdditionalHeader()

void ParallelIoPNetcdf::b_addAdditionalHeader ( )
private
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2013-04-12

Definition at line 561 of file parallelio_pnetcdf.cpp.

561 {
562 TRACE();
563
564 b_ncRedef();
565
566 // For a newly created or modified file, add some meta information
567 if(m_fileMode != PIO_READ) {
568 const MInt maxNoChars = 256;
569
570 // Get all meta-data on root process & communicate, since some of the C
571 // functions might not be thread-safe (at least getpwuid() is not)
572 MChar user[maxNoChars];
573 MChar host[maxNoChars];
574 MChar dir[maxNoChars];
575 MChar exec[maxNoChars];
576 MChar date[maxNoChars];
577
578 // Create object & initialize data
579 fill(user, user + maxNoChars, '\0');
580 fill(host, host + maxNoChars, '\0');
581 fill(dir, dir + maxNoChars, '\0');
582 fill(exec, exec + maxNoChars, '\0');
583 fill(date, date + maxNoChars, '\0');
584
585 if(m_domainId == 0) {
586 // Gets the current username
587 passwd* p;
588 p = getpwuid(getuid());
589 if(p) {
590 strncpy(user, p->pw_name, maxNoChars - 1);
591 } else {
592 strncpy(user, "n/a", maxNoChars - 1);
593 }
594
595 // Gets the current hostname
596 gethostname(host, maxNoChars - 1);
597 strcpy(&host[strlen(host)], " (");
598 strcpy(&host[strlen(host)], XSTRINGIFY(MAIA_HOST_STRING));
599 strcpy(&host[strlen(host)], ")");
600
601// Gets the current directory
602#ifdef MAIA_GCC_COMPILER
603#pragma GCC diagnostic push
604#pragma GCC diagnostic ignored "-Wunused-result"
605#endif
606
607 getcwd(dir, maxNoChars - 1);
608
609 readlink("/proc/self/exe", exec, maxNoChars - 1);
610
611#ifdef MAIA_GCC_COMPILER
612#pragma GCC diagnostic pop
613#endif
614
615 // Get the current time and write it to rawTime
616 time_t rawTime;
617 time(&rawTime);
618
619 // Convert to time struct
620 tm* timeInfo;
621 timeInfo = localtime(&rawTime);
622
623 // Format time to string and save to buffer
624 strftime(date, maxNoChars, "%Y-%m-%d %H:%M:%S", timeInfo);
625 }
626
627 // Pack data
628 const MInt noItems = 5;
629 MChar buffer[noItems * maxNoChars];
630 memcpy(buffer + 0 * maxNoChars, user, maxNoChars);
631 memcpy(buffer + 1 * maxNoChars, host, maxNoChars);
632 memcpy(buffer + 2 * maxNoChars, dir, maxNoChars);
633 memcpy(buffer + 3 * maxNoChars, exec, maxNoChars);
634 memcpy(buffer + 4 * maxNoChars, date, maxNoChars);
635
636 // Broadcast time from rank 0 to ensure that every rank has the same
637 // information
638 MPI_Bcast(&buffer, noItems * maxNoChars, MPI_CHAR, 0, m_mpiComm, AT_, "buffer");
639
640 // Unpack data
641 memcpy(user, buffer + 0 * maxNoChars, maxNoChars);
642 memcpy(host, buffer + 1 * maxNoChars, maxNoChars);
643 memcpy(dir, buffer + 2 * maxNoChars, maxNoChars);
644 memcpy(exec, buffer + 3 * maxNoChars, maxNoChars);
645 memcpy(date, buffer + 4 * maxNoChars, maxNoChars);
646 MString version = MString(XSTRINGIFY(MAIA_VERSION_STRING));
647 MString build = MString(XSTRINGIFY(MAIA_COMPILER_STRING)) + " " + MString(XSTRINGIFY(MAIA_BUILD_TYPE_STRING)) + " ("
648 + MString(XSTRINGIFY(MAIA_COMPILER_VERSION_STRING)) + ")";
649
651 // Add file attributes only needed for creation
652 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(MString(user), "_meta_creation_user");
653 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(MString(host), "_meta_creation_host");
654 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(MString(dir), "_meta_creation_directory");
655 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(MString(exec), "_meta_creation_executable");
657 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(MString(date), "_meta_creation_date");
658 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(version, "_meta_creation_revision");
659 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(build, "_meta_creation_build");
660 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(ncmpi_inq_libvers(), "_meta_creation_pnetcdf_version");
661 } else if(m_fileMode == PIO_APPEND) {
662 // Add file attributes that should be set at each modification
663 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(MString(user), "_meta_lastModified_user");
664 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(MString(host), "_meta_lastModified_host");
665 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(MString(dir), "_meta_lastModified_directory");
666 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(MString(exec), "_meta_lastModified_executable");
667 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(m_noDomains, "_meta_lastModified_noDomains");
668 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(MString(date), "_meta_lastModified_date");
669 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(version, "_meta_lastModified_revision");
670 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(build, "_meta_lastModified_build");
671 ParallelIoBase<ParallelIoPNetcdf>::setAttribute(ncmpi_inq_libvers(), "_meta_lastModified_pnetcdf_version");
672 }
673 }
674}
void setAttribute(const T &value, const MString &name, const MString &datasetName="")
Set a file or dataset attribute. [MPI]
Definition: parallelio.h:1438
std::basic_string< char > MString
Definition: maiatypes.h:55
char MChar
Definition: maiatypes.h:56
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
constexpr std::underlying_type< FcCell >::type p(const FcCell property)
Converts property name to underlying integer value.

◆ b_defineArray() [1/2]

void ParallelIoPNetcdf::b_defineArray ( maiabd_type  type,
const MString name,
size_type  noDims,
size_type totalCount 
)
private
Author
Ramandeep Jain (HiWi) raman.nosp@m.deep.nosp@m.jain@.nosp@m.gmai.nosp@m.l.com, Konstantin Froehlich
Date
2015-Mar
Parameters
[in]typeData type of the array (may be either maia::parallel_io::PIO_FLOAT or maia::parallel_io::PIO_INT or maia::parallel_io::PIO_STRING or maia::parallel_io::PIO_UCHAR).
[in]nameName of the array (must not be empty).
[in]noDimsNumber of array dimensions
[in]totalCountTotal size of the array in each dimension.

Definition at line 717 of file parallelio_pnetcdf.cpp.

717 {
718 TRACE();
719
720 b_ncRedef();
721
722 // Fill values
723 MFloat fillValueFloat = MFloatNaN;
724 MInt fillValueInt = std::numeric_limits<MInt>::max();
725 MLong fillValueLong = std::numeric_limits<MLong>::max();
726 MUlong fillValueUlong = std::numeric_limits<MUlong>::max();
727 MUchar fillValueChar = std::numeric_limits<MChar>::max();
728 MUchar fillValueUchar = std::numeric_limits<MUchar>::max();
729
730 // Determine NC data type
731 nc_type dataType;
732 [[maybe_unused]] void* fillValue = nullptr;
733 if(type == PIO_FLOAT) {
734 dataType = NC_DOUBLE;
735 fillValue = &fillValueFloat;
736 } else if(type == PIO_INT) {
737 dataType = NC_INT;
738 fillValue = &fillValueInt;
739 } else if(type == PIO_LONG) {
740 dataType = NC_INT64;
741 fillValue = &fillValueLong;
742 } else if(type == PIO_STRING) {
743 dataType = NC_CHAR;
744 fillValue = &fillValueChar;
745 } else if(type == PIO_UCHAR) {
746 dataType = NC_UBYTE;
747 fillValue = &fillValueUchar;
748 } else if(type == PIO_ULONGLONG) {
749 dataType = NC_UINT64;
750 fillValue = &fillValueUlong;
751 } else {
752 TERMM(1, "Invalid ParallelIo data type!");
753 }
754
755 // Determine whether (a) new dimension(s) needs to be created and if yes,
756 // create it
757 MInt dimId, status, varId;
758 MIntScratchSpace dimIds(noDims, FUN_, "dimIds");
759 for(size_type dId = 0; dId < noDims; dId++) {
760 if(b_ncDimensions.count(totalCount[dId]) != 0u) {
761 // Dimension was found, now get its id
762 dimId = b_ncDimensions.find(totalCount[dId])->second.dimId;
763 } else {
764 // Get next dimension name
765 MInt maxUsedDimensionNo = -1;
766 for(NcDimMap::const_iterator it = b_ncDimensions.begin(); it != b_ncDimensions.end(); it++) {
767 maxUsedDimensionNo = max(maxUsedDimensionNo, it->second.dimNo);
768 }
769
770 // Set dimension name according to found dimension numbers
771 const MInt dimNo = maxUsedDimensionNo + 1;
772 const MString dimensionName = "dim" + to_string(dimNo);
773
774 // Define new dimension in the file and store its information
775 status = ncmpi_def_dim(b_ncId, dimensionName.c_str(), totalCount[dId], &dimId);
776 b_error(status, dimensionName, AT_);
777
778 // Add dimension to map
779 b_ncDimensions[totalCount[dId]] = NcDimProxy();
780 b_ncDimensions[totalCount[dId]].dimNo = dimNo;
781 b_ncDimensions[totalCount[dId]].dimId = dimId;
782 }
783 dimIds[dId] = dimId;
784 }
785
786 // Define dimension and variable in the file
787 status = ncmpi_def_var(b_ncId, name.c_str(), dataType, noDims, &dimIds[0], &varId);
788 b_error(status, name, AT_);
789
790#if MAIA_NCMPI_FILL_VARIABLES
791 // NOTE: set fill mode and default fill value to detect unwritten variables/detect IO errors
792 status = ncmpi_def_var_fill(b_ncId, varId, 0, fillValue);
793 b_error(status, name, AT_);
794#endif
795}
MLong size_type
Type used for all size- and offset-related values.
Definition: parallelio.h:123
This class is a ScratchSpace.
Definition: scratch.h:758
unsigned char MUchar
Definition: maiatypes.h:57
double MFloat
Definition: maiatypes.h:52
int64_t MLong
Definition: maiatypes.h:64
uint64_t MUlong
Definition: maiatypes.h:65
const MInt PIO_STRING
Definition: parallelio.h:52
const MInt PIO_ULONGLONG
Definition: parallelio.h:56
const MInt PIO_UCHAR
Definition: parallelio.h:54
const MInt PIO_INT
Definition: parallelio.h:48
const MInt PIO_LONG
Definition: parallelio.h:50
const MInt PIO_FLOAT
Definition: parallelio.h:46

◆ b_defineArray() [2/2]

void ParallelIoPNetcdf::b_defineArray ( maiabd_type  ,
const MString ,
const MString ,
const  size_type,
const size_type  
)
private

◆ b_defineScalar()

void ParallelIoPNetcdf::b_defineScalar ( maiabd_type  type,
const MString name 
)
private
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]typeData type of the scalar (may be either maia::parallel_io::PIO_FLOAT or maia::parallel_io::PIO_INT).
[in]nameName of the scalar (must not be empty).

Definition at line 808 of file parallelio_pnetcdf.cpp.

808 {
809 TRACE();
810
811 // Choose the HDF5 specific native datatype or a string of variable length
812 b_ncRedef();
813
814 // Fill values
815 MFloat fillValueFloat = MFloatNaN;
816 MInt fillValueInt = std::numeric_limits<MInt>::max();
817 MLong fillValueLong = std::numeric_limits<MLong>::max();
818 MUlong fillValueUlong = std::numeric_limits<MUlong>::max();
819 MUchar fillValueChar = std::numeric_limits<MChar>::max();
820 MUchar fillValueUchar = std::numeric_limits<MUchar>::max();
821
822 // Determine NC data type
823 nc_type dataType;
824 [[maybe_unused]] void* fillValue = nullptr;
825 if(type == PIO_FLOAT) {
826 dataType = NC_DOUBLE;
827 fillValue = &fillValueFloat;
828 } else if(type == PIO_INT) {
829 dataType = NC_INT;
830 fillValue = &fillValueInt;
831 } else if(type == PIO_LONG) {
832 dataType = NC_INT64;
833 fillValue = &fillValueLong;
834 } else if(type == PIO_STRING) {
835 dataType = NC_CHAR;
836 fillValue = &fillValueChar;
837 } else if(type == PIO_UCHAR) {
838 dataType = NC_UBYTE;
839 fillValue = &fillValueUchar;
840 } else if(type == PIO_ULONGLONG) {
841 dataType = NC_UINT64;
842 fillValue = &fillValueUlong;
843 } else {
844 TERMM(1, "Invalid ParallelIo data type!");
845 }
846
847 // For a scalar, only the variable needs to be defined
848 MInt varId;
849 MInt status = ncmpi_def_var(b_ncId, name.c_str(), dataType, 0, nullptr, &varId);
850 b_error(status, name, AT_);
851
852#if MAIA_NCMPI_FILL_VARIABLES
853 // NOTE: set fill mode and default fill value to detect unwritten variables/detect IO errors
854 status = ncmpi_def_var_fill(b_ncId, varId, 0, fillValue);
855 b_error(status, name, AT_);
856#endif
857}

◆ b_error()

void ParallelIoPNetcdf::b_error ( MInt  status,
const MString name,
const MString location 
)
staticprivate
Author
Ramandeep Jain (HiWi) raman.nosp@m.deep.nosp@m.jain@.nosp@m.gmai.nosp@m.l.com, Konstantin Froehlich
Date
2015-Mar
Parameters
[in]statusHDF5 status code (negative if unexspected behaviour)
[in]nameName of the file / variable / attribute
[in]locationFunction which called the error message
Returns
The status message if there was an error, or zero if everything is ok.

Definition at line 2100 of file parallelio_pnetcdf.cpp.

2100 {
2101 // TRACE(); //<- this function is called to often and not really important.
2102 if(status != NC_NOERR) {
2103 cerr << endl;
2104 cerr << "*** ERROR in parallelio_pnetcdf ***" << endl;
2105 cerr << "NetCDF error in '" << location << "'" << endl;
2106 cerr << "NetCDF returns status " << status << ": " << ncmpi_strerror(status) << endl;
2107 cerr << "The file/variable/attribute in question was: " << name << endl;
2108 cerr << endl;
2109 TERMM(1, "NetCDF error in ParallelIo.");
2110 }
2111}
const MString & location
Definition: functions.h:37

◆ b_fileExt()

const MString ParallelIoPNetcdf::b_fileExt ( )
staticprivate
Author
Konstantin Froehlich
Date
2015-10-21

Definition at line 349 of file parallelio_pnetcdf.cpp.

349 {
350 TRACE();
351
352 return ".Netcdf";
353}

◆ b_getAttribute() [1/2]

template<>
void ParallelIoPNetcdf::b_getAttribute ( MString *const  value,
const MString name,
const MString datasetName,
const size_type  totalCount 
)
Author
Ramandeep Jain (HiWi) raman.nosp@m.deep.nosp@m.jain@.nosp@m.gmai.nosp@m.l.com, Konstantin Froehlich
Date
2015-Mar
Parameters
[out]valueAttribute value
[in]nameAttribute name
[in]datasetNameDataset (scalar or array) to which the attribute is attached (leave empty for file attributes).

Definition at line 2030 of file parallelio_pnetcdf.cpp.

2031 {
2032 TRACE();
2033
2034 if(totalCount > 1) {
2035 mTerm(1, AT_, "Array of strings attributes not yet supported.");
2036 }
2037
2038 // Determine variable id
2039 MInt varId;
2040 if(datasetName.empty()) {
2041 varId = NC_GLOBAL;
2042 } else {
2043 MInt status = ncmpi_inq_varid(b_ncId, datasetName.c_str(), &varId);
2044 b_error(status, datasetName, AT_);
2045 }
2046
2047 // Get attribute length
2048 MInt status;
2049 MPI_Offset length;
2050 status = ncmpi_inq_attlen(b_ncId, varId, name.c_str(), &length);
2051 b_error(status, datasetName + "::" + name, AT_);
2052
2053 // Read attribute
2054 ScratchSpace<MChar> tmpScratch(length, FUN_, "tmpScratch");
2055 status = ncmpi_get_att_text(b_ncId, varId, name.c_str(), tmpScratch.data());
2056 b_error(status, datasetName + "::" + name, AT_);
2057 value->assign(tmpScratch.data(), length);
2058}

◆ b_getAttribute() [2/2]

template<class T >
void ParallelIoPNetcdf::b_getAttribute ( T *const  value,
const MString name,
const MString datasetName = "",
const size_type  totalCount = 1 
)
private
Author
Ramandeep Jain (HiWi) raman.nosp@m.deep.nosp@m.jain@.nosp@m.gmai.nosp@m.l.com, Konstantin Froehlich
Date
2015-Mar
Parameters
[out]valueAttribute value
[in]nameAttribute name
[in]datasetNameDataset (scalar or array) to which the attribute is attached (leave empty for file attributes).

Definition at line 1987 of file parallelio_pnetcdf.cpp.

1988 {
1989 TRACE();
1990
1991 // Determine variable id
1992 MInt varId;
1993 if(datasetName.empty()) {
1994 varId = NC_GLOBAL;
1995 } else {
1996 MInt status = ncmpi_inq_varid(b_ncId, datasetName.c_str(), &varId);
1997 b_error(status, datasetName, AT_);
1998 }
1999
2000 // Get attribute length
2001 MInt status;
2002 MPI_Offset length;
2003 status = ncmpi_inq_attlen(b_ncId, varId, name.c_str(), &length);
2004 b_error(status, datasetName + "::" + name, AT_);
2005
2006 if(length != (MPI_Offset)totalCount) {
2007 TERMM(1, "Requested attribute (" + name + ") has different number of elements (" + to_string(length)
2008 + ") than specified (" + to_string(totalCount)
2009 + "). Use getAttributeCount() to query number of elements first");
2010 }
2011
2012 // Read attribute
2013 status = pnetcdf_traits<T>::ncmpi_get_att_type(b_ncId, varId, name.c_str(), value);
2014 b_error(status, datasetName + "::" + name, AT_);
2015}

◆ b_getAttributeCount()

void ParallelIoPNetcdf::b_getAttributeCount ( const MString name,
size_type totalCount,
const MString datasetName = "" 
)
private

Definition at line 2061 of file parallelio_pnetcdf.cpp.

2061 {
2062 TRACE();
2063
2064 // Determine variable id
2065 MInt varId;
2066 if(datasetName.empty()) {
2067 varId = NC_GLOBAL;
2068 } else {
2069 MInt status = ncmpi_inq_varid(b_ncId, datasetName.c_str(), &varId);
2070 b_error(status, datasetName, AT_);
2071 }
2072
2073 // Get attribute length
2074 MInt status;
2075 MPI_Offset length;
2076 status = ncmpi_inq_attlen(b_ncId, varId, name.c_str(), &length);
2077 b_error(status, datasetName + "::" + name, AT_);
2078
2079 *totalCount = (size_type)length;
2080}

◆ b_getAttributeType()

MInt ParallelIoPNetcdf::b_getAttributeType ( const MString name,
const MString datasetName = "" 
)
private
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.

Definition at line 1174 of file parallelio_pnetcdf.cpp.

1174 {
1175 TRACE();
1176
1177 // Determine variable id
1178 MInt varId;
1179 if(datasetName.empty()) {
1180 varId = NC_GLOBAL;
1181 } else {
1182 MInt status = ncmpi_inq_varid(b_ncId, datasetName.c_str(), &varId);
1183 b_error(status, datasetName, AT_);
1184 }
1185
1186 // Get attribute type
1187 nc_type ncType;
1188 MInt status = ncmpi_inq_atttype(b_ncId, varId, name.c_str(), &ncType);
1189 b_error(status, name, AT_);
1190
1191 // Translate NC type to ParallelIo type
1192 MInt typeId;
1193 if(ncType == NC_INT) {
1194 typeId = PIO_INT;
1195 } else if(ncType == NC_DOUBLE) {
1196 typeId = PIO_FLOAT;
1197 } else if(ncType == NC_INT64) {
1198 typeId = PIO_LONG;
1199 } else if(ncType == NC_CHAR) {
1200 typeId = PIO_STRING;
1201 } else if(ncType == NC_UBYTE) {
1202 typeId = PIO_UCHAR;
1203 } else if(ncType == NC_UINT64) {
1204 typeId = PIO_ULONGLONG;
1205 } else {
1206 typeId = PIO_UNKNOWN_TYPE;
1207 }
1208
1209 return typeId;
1210}
const MInt PIO_UNKNOWN_TYPE
Definition: parallelio.h:44

◆ b_getDatasetNames() [1/2]

void ParallelIoPNetcdf::b_getDatasetNames ( std::vector< MString > &  names,
const MString path 
)
private

Definition at line 925 of file parallelio_pnetcdf.cpp.

925 {
926 TRACE();
927
928 (void)path;
929 (void)names;
930 mTerm(1, AT_, "Group functionality not supported by PNetcdf backend");
931}

◆ b_getDatasetNames() [2/2]

void ParallelIoPNetcdf::b_getDatasetNames ( std::vector< MString > &  names,
const size_type  dimension 
)
private
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]dimensionMatch only datasets with given number of dimensions, default -1 which will return all datasets, 0 will only return scalars, 1 arrays with one dimension, ...
[out]namesList of dataset names.

Definition at line 1020 of file parallelio_pnetcdf.cpp.

1020 {
1021 TRACE();
1022
1023 // Erase vector contents
1024 vector<MString>().swap(names);
1025
1026 MInt status;
1027
1028 // Get number of variable in file
1029 MInt noVars;
1030 status = ncmpi_inq_nvars(b_ncId, &noVars);
1031 b_error(status, m_fileName, AT_);
1032
1033 // Check each variable if
1034 // - has the specified number of dimensions or return all if dimension=-1
1035 // (default)
1036 for(MInt i = 0; i < noVars; i++) {
1037 MInt noDimensions;
1038 status = ncmpi_inq_varndims(b_ncId, i, &noDimensions);
1039 b_error(status, m_fileName, AT_);
1040
1041 if(noDimensions == dimension || dimension == -1) {
1042 // If this is a data file, get data file-specific name, otherwise just get
1043 // the variable name
1044 MChar varname[NC_MAX_NAME + 1];
1045 status = ncmpi_inq_varname(b_ncId, i, varname);
1046 b_error(status, m_fileName, AT_);
1047 names.emplace_back(varname);
1048 }
1049 }
1050}

◆ b_getDatasetNoDims() [1/2]

ParallelIo::size_type ParallelIoPNetcdf::b_getDatasetNoDims ( const MString name)
private
Author
Konstantin Froehlich
Date
2015-10-09
Parameters
[in]nameDataset name for which the number of dimensions is returned

The number of dimensions of a dataset is returned, 0 is a scalar, 1 an array with one dimension, ...

Returns
Number of dimensions of the specified dataset

Definition at line 1066 of file parallelio_pnetcdf.cpp.

1066 {
1067 TRACE();
1068
1069 // Get variable id
1070 MInt varId;
1071 MInt status = ncmpi_inq_varid(b_ncId, name.c_str(), &varId);
1072 b_error(status, name, AT_);
1073
1074 // Get number of variable dimensions
1075 MInt noDims;
1076 status = ncmpi_inq_varndims(b_ncId, varId, &noDims);
1077 b_error(status, name, AT_);
1078
1079 return noDims;
1080}

◆ b_getDatasetNoDims() [2/2]

ParallelIo::size_type ParallelIoPNetcdf::b_getDatasetNoDims ( const MString path,
const MString name 
)
private

Definition at line 941 of file parallelio_pnetcdf.cpp.

941 {
942 TRACE();
943
944 (void)path;
945 (void)name;
946 mTerm(1, AT_, "Group functionality not supported by PNetcdf backend");
947 return -1;
948}

◆ b_getDatasetSize() [1/2]

void ParallelIoPNetcdf::b_getDatasetSize ( const MString name,
const MString path,
size_type  noDims,
size_type data 
)
private

Definition at line 950 of file parallelio_pnetcdf.cpp.

950 {
951 TRACE();
952
953 (void)path;
954 (void)name;
955 (void)noDims;
956 (void)data;
957 mTerm(1, AT_, "Group functionality not supported by PNetcdf backend");
958}

◆ b_getDatasetSize() [2/2]

ParallelIo::size_type ParallelIoPNetcdf::b_getDatasetSize ( const MString name,
const size_type  dimensionId = 0 
)
private
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.
Returns
Total number of elements of the array in the given dimension.

Definition at line 1094 of file parallelio_pnetcdf.cpp.

1094 {
1095 TRACE();
1096
1097 // Get variable id
1098 MInt varId;
1099 MInt status = ncmpi_inq_varid(b_ncId, name.c_str(), &varId);
1100 b_error(status, name, AT_);
1101
1102 // Get number of array dimensions
1103 size_type noDims = b_getDatasetNoDims(name);
1104
1105 // Get variable dimension
1106 MIntScratchSpace dimId(noDims, FUN_, "dimId");
1107 status = ncmpi_inq_vardimid(b_ncId, varId, &dimId[0]);
1108 b_error(status, name, AT_);
1109
1110 // Get variable size
1111 MPI_Offset arraySize;
1112 status = ncmpi_inq_dimlen(b_ncId, dimId[dimensionId], &arraySize);
1113 b_error(status, name, AT_);
1114
1115 return static_cast<size_type>(arraySize);
1116}
size_type b_getDatasetNoDims(const MString &name)
Get number of dimensions of a dataset with the given name.

◆ b_getDatasetType()

MInt ParallelIoPNetcdf::b_getDatasetType ( const MString name)
private
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 dataset to check.
Returns
The data type as defined in the namespace maia::parallel_io.

Definition at line 972 of file parallelio_pnetcdf.cpp.

972 {
973 TRACE();
974
975 // Get variable id
976 MInt varId;
977 MInt status = ncmpi_inq_varid(b_ncId, name.c_str(), &varId);
978 b_error(status, name, AT_);
979
980 // Get variable type
981 nc_type ncType;
982 status = ncmpi_inq_vartype(b_ncId, varId, &ncType);
983 b_error(status, name, AT_);
984
985 // Translate NC type to ParallelIo type
986 MInt typeId;
987 if(ncType == NC_INT) {
988 typeId = PIO_INT;
989 } else if(ncType == NC_DOUBLE) {
990 typeId = PIO_FLOAT;
991 } else if(ncType == NC_INT64) {
992 typeId = PIO_LONG;
993 } else if(ncType == NC_CHAR) {
994 typeId = PIO_STRING;
995 } else if(ncType == NC_UBYTE) {
996 typeId = PIO_UCHAR;
997 } else if(ncType == NC_UINT64) {
998 typeId = PIO_ULONGLONG;
999 } else {
1000 typeId = PIO_UNKNOWN_TYPE;
1001 }
1002
1003 return typeId;
1004}

◆ b_getGroupNames()

void ParallelIoPNetcdf::b_getGroupNames ( std::vector< MString > &  names,
const MString path 
)
private

Definition at line 933 of file parallelio_pnetcdf.cpp.

933 {
934 TRACE();
935
936 (void)path;
937 (void)names;
938 mTerm(1, AT_, "Group functionality not supported by PNetcdf backend");
939}

◆ b_hasAttribute()

MBool ParallelIoPNetcdf::b_hasAttribute ( const MString name,
const MString datasetName = "" 
)
private
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 1131 of file parallelio_pnetcdf.cpp.

1131 {
1132 TRACE();
1133
1134 // Determine variable id
1135 MInt varId;
1136 if(datasetName.empty()) {
1137 varId = NC_GLOBAL;
1138 } else {
1139 MInt status = ncmpi_inq_varid(b_ncId, datasetName.c_str(), &varId);
1140 b_error(status, datasetName, AT_);
1141 }
1142
1143 // Determine if attribute exists
1144 MInt attId;
1145 MInt status = ncmpi_inq_attid(b_ncId, varId, name.c_str(), &attId);
1146
1147 // Generate return falue
1148 MBool attributeExists;
1149 if(status == NC_NOERR) {
1150 attributeExists = true;
1151 } else if(status == NC_ENOTATT) {
1152 attributeExists = false;
1153 } else {
1154 attributeExists = false;
1155 b_error(status, name, AT_);
1156 }
1157
1158 return attributeExists;
1159}
bool MBool
Definition: maiatypes.h:58

◆ b_hasDataset() [1/2]

MBool ParallelIoPNetcdf::b_hasDataset ( const MString name,
const MString path 
)
private

Definition at line 916 of file parallelio_pnetcdf.cpp.

916 {
917 TRACE();
918
919 (void)path;
920 (void)name;
921 mTerm(1, AT_, "Group functionality not supported by PNetcdf backend");
922 return false;
923}

◆ b_hasDataset() [2/2]

MBool ParallelIoPNetcdf::b_hasDataset ( const MString name,
const size_type  noDimensions 
)
private
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]noDimensionsNumber of dimensions of the dataset to match, default -1 which is any dimensionality, 0 will only check scalars, 1 arrays with one dimension, ...
Returns
True if an dataset of the given name with the given number of dimensions exists.

Definition at line 878 of file parallelio_pnetcdf.cpp.

878 {
879 TRACE();
880
881 // Get number of variable in file
882 MInt noVars;
883 MInt status = ncmpi_inq_nvars(b_ncId, &noVars);
884 b_error(status, m_fileName, AT_);
885
886 // Check each variable if
887 // - names match
888 // - has the given dimension (or match any dimension if -1 which is default)
889 MBool varExists = false;
890 for(MInt i = 0; i < noVars; i++) {
891 MChar varname[NC_MAX_NAME + 1];
892 status = ncmpi_inq_varname(b_ncId, i, varname);
893 b_error(status, m_fileName, AT_);
894
895 MInt nDims;
896 status = ncmpi_inq_varndims(b_ncId, i, &nDims);
897 b_error(status, varname, AT_);
898
899 if(name == varname && (nDims == noDimensions || noDimensions == -1)) {
900 varExists = true;
901 break;
902 }
903 }
904
905 return varExists;
906}

◆ b_hasObject()

MBool ParallelIoPNetcdf::b_hasObject ( const MString path)
private

Definition at line 908 of file parallelio_pnetcdf.cpp.

908 {
909 TRACE();
910
911 // pnetcdf doesn't support groups (hdf5 feature), therefore this
912 // function can only check for datasets
913 return b_hasDataset(path, -1);
914}
MBool b_hasDataset(const MString &name, const size_type dimension)
Check if the file contains a dataset with the given name and number of dimensions.

◆ b_isValidFile()

MBool ParallelIoPNetcdf::b_isValidFile ( const MString name,
const MPI_Comm &  mpiComm 
)
staticprivate
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 file that is to be checked.
[in]mpiCommMPI communicator to be used.

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

Definition at line 322 of file parallelio_pnetcdf.cpp.

322 {
323 TRACE();
324
325 MBool returnValue;
326
327 MInt status, fileId;
328 status = ncmpi_open(mpiComm, name.c_str(), NC_NOWRITE, globalMpiInfo(), &fileId);
329
330 if(status == NC_NOERR) {
331 returnValue = true;
332 status = ncmpi_close(fileId);
333 b_error(status, name, AT_);
334 } else {
335 returnValue = false;
336 }
337
338 return returnValue;
339}

◆ b_ncEndDef()

void ParallelIoPNetcdf::b_ncEndDef ( )
private
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

Definition at line 502 of file parallelio_pnetcdf.cpp.

502 {
503 TRACE();
504
505 if(!b_ncDataMode) {
506 MInt status = ncmpi_enddef(b_ncId);
507 b_error(status, m_fileName, AT_);
508 b_ncDataMode = true;
509 }
510}

◆ b_ncRedef()

void ParallelIoPNetcdf::b_ncRedef ( )
private
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

Definition at line 519 of file parallelio_pnetcdf.cpp.

519 {
520 TRACE();
521
522 if(b_ncDataMode) {
523 MInt status = ncmpi_redef(b_ncId);
524 b_error(status, m_fileName, AT_);
525 b_ncDataMode = false;
526 }
527}

◆ b_printFileHints()

void ParallelIoPNetcdf::b_printFileHints ( )
private

Source: https://trac.mcs.anl.gov/projects/parallel-netcdf/browser/trunk/examples/C/hints.c

Definition at line 1218 of file parallelio_pnetcdf.cpp.

1218 {
1219 TRACE();
1220
1221 char value[MPI_MAX_INFO_VAL];
1222 MInt status, len, flag;
1223 MPI_Offset header_size, header_extent;
1224 MPI_Offset h_align = -1, v_align = -1, h_chunk = -1;
1225 MPI_Info info_used;
1226
1227 // Get header size
1228 status = ncmpi_inq_header_size(b_ncId, &header_size);
1229 b_error(status, m_fileName, AT_);
1230
1231 // Get maximum size of header (without the need to move data)
1232 status = ncmpi_inq_header_extent(b_ncId, &header_extent);
1233 b_error(status, m_fileName, AT_);
1234
1235 // Get file MPI information
1236 status = ncmpi_inq_file_info(b_ncId, &info_used);
1237 b_error(status, m_fileName, AT_);
1238
1239 // Get header align size (in bytes)
1240 MPI_Info_get_valuelen(info_used, "nc_header_align_size", &len, &flag, AT_);
1241 if(flag != 0) {
1242 MPI_Info_get(info_used, "nc_header_align_size", len + 1, value, &flag, AT_);
1243 h_align = strtoll(value, nullptr, 10);
1244 }
1245
1246 // Get variable align size (in bytes)
1247 MPI_Info_get_valuelen(info_used, "nc_var_align_size", &len, &flag, AT_);
1248 if(flag != 0) {
1249 MPI_Info_get(info_used, "nc_var_align_size", len + 1, value, &flag, AT_);
1250 v_align = strtoll(value, nullptr, 10);
1251 }
1252
1253 // Get header read chuck size (in bytes)
1254 MPI_Info_get_valuelen(info_used, "nc_header_read_chunk_size", &len, &flag, AT_);
1255 if(flag != 0) {
1256 MPI_Info_get(info_used, "nc_header_read_chunk_size", len + 1, value, &flag, AT_);
1257 h_chunk = strtoll(value, nullptr, 10);
1258 }
1259
1260 MPI_Info_free(&info_used, AT_);
1261
1262 // Output file hint information
1263 std::cerr << "##### PNetCDF file hints #####" << std::endl;
1264
1265 if(h_align == -1) {
1266 std::cerr << "nc_header_align_size is NOT set" << std::endl;
1267 } else {
1268 std::cerr << "nc_header_align_size set to = " << h_align << std::endl;
1269 }
1270
1271 if(v_align == -1) {
1272 std::cerr << "nc_var_align_size is NOT set" << std::endl;
1273 } else {
1274 std::cerr << "nc_var_align_size set to = " << v_align << std::endl;
1275 }
1276 if(h_chunk == -1) {
1277 std::cerr << "nc_header_read_chunk_size is NOT set" << std::endl;
1278 } else {
1279 std::cerr << "nc_header_read_chunk_size set to = " << h_chunk << std::endl;
1280 }
1281
1282 std::cerr << "header size = " << header_size << std::endl;
1283 std::cerr << "header extent = " << header_extent << std::endl;
1284
1285 // Check free header space in append mode
1286 if(m_fileMode == PIO_APPEND) {
1287 const MPI_Offset header_free = header_extent - header_size;
1288 std::cerr << "header free space (append mode) = " << header_free << std::endl;
1289 if(header_free < 1024) {
1290 std::cerr << "WARNING: ParallelIoPNetcdf file in append mode has less than 1KB of free header "
1291 "space. If this space is used up by adding new data to the header the entire data "
1292 "file has to be moved to make room for the new definitions which may cause MPI I/O "
1293 "errors."
1294 << std::endl;
1295 }
1296 }
1297
1298 std::cerr << "##### PNetCDF file hints #####" << std::endl;
1299}
int MPI_Info_get(MPI_Info info, const char *key, int valuelen, char *value, int *flag, const MString &name)
same as MPI_Info_get
int MPI_Info_free(MPI_Info *info, const MString &name)
same as MPI_Info_free
int MPI_Info_get_valuelen(MPI_Info info, const char *key, int *valuelen, int *flag, const MString &name)
same as MPI_Info_get_valuelen

◆ b_readArray() [1/3]

template<>
void ParallelIoPNetcdf::b_readArray ( MString array,
const MString name,
const size_type  noDims,
MPI_Offset *  start,
MPI_Offset *  count,
MPI_Offset  memoryStride,
const size_type  noChunks,
MPI_Offset  diskStride 
)
Author
Ansgar Niemoeller, Konstantin Froehlich
Date
2015-08-19
Parameters
[in]arrayPointer to data in memory.
[in]nameName of the dataset.
[in]noDimsNumber of array dimensions
[in]startStart offset in the file (i.e. the domain should start reading from here).
[in]countNumber of elements read from this domain.
[in]memoryStrideStride in memory between data values.

This method is the exact counterpart to h5WriteArray().

Definition at line 1726 of file parallelio_pnetcdf.cpp.

1728 {
1729 TRACE();
1730
1731 b_ncEndDef();
1732
1733 // Get variable id
1734 MInt varId;
1735 MInt status = ncmpi_inq_varid(b_ncId, name.c_str(), &varId);
1736 b_error(status, name, AT_);
1737
1738 // Determine total data count
1739 size_type totalCount = 1;
1740 for(size_type d = 0; d < noDims - 1; d++) {
1741 totalCount *= count[d];
1742 }
1743
1744 // Determine length of one string
1745 size_type strLen = count[noDims - 1];
1746
1747 // Create temporary storage space if needed and set data pointers
1748 MInt tmpScratchSize = (memoryStride == 1 ? 1 : totalCount);
1749 ScratchSpace<MString> tmpScratch(tmpScratchSize, FUN_, "tmpStorage");
1750 MString* data = nullptr;
1751 if(memoryStride == 1) {
1752 data = array;
1753 } else {
1754 data = tmpScratch.data();
1755 }
1756
1757 // Create buffer for reading
1758 size_type nCount = totalCount * strLen;
1759 ScratchSpace<char> buf(nCount, FUN_, "buf");
1760
1761 // Read array
1762 if(noChunks == 1) {
1763 // If number of chunks is one, read everything at once
1764 if(diskStride == 1) {
1765 status = ncmpi_get_vara_text_all(b_ncId, varId, start, count, &buf[0]);
1766 } else {
1767 status = ncmpi_get_vars_text_all(b_ncId, varId, start, count, &diskStride, &buf[0]);
1768 }
1769 b_error(status, name, AT_);
1770
1771 // Extract strings from buffer
1772 for(size_type i = 0; i < totalCount; i++) {
1773 MString tmp;
1774 tmp.append(&buf[i * strLen], strLen);
1775 data[i].append(tmp.c_str(), 0, strLen);
1776 }
1777 } else {
1778 // Read in chunks
1779 MPI_Offset chunkSize = count[0] / noChunks;
1780 if(count[0] % noChunks > 0) {
1781 chunkSize += 1;
1782 }
1783
1784 ScratchSpace<MPI_Offset> start_(noDims, FUN_, "start_");
1785 ScratchSpace<MPI_Offset> count_(noDims, FUN_, "count_");
1786
1787 std::copy(start, start + noDims, &start_[0]);
1788 std::copy(count, count + noDims, &count_[0]);
1789
1790 for(size_type i = 0; i < noChunks; i++) {
1791 start_[0] = min(start[0] + i * chunkSize * diskStride, count[0] * diskStride - 1);
1792 count_[0] = max(min(chunkSize, count[0] - i * chunkSize), 0ll);
1793 if(diskStride == 1) {
1794 status = ncmpi_get_vara_text_all(b_ncId, varId, &start_[0], &count_[0], &buf[0]);
1795 } else {
1796 status = ncmpi_get_vars_text_all(b_ncId, varId, &start_[0], &count_[0], &diskStride, &buf[0]);
1797 }
1798
1799 b_error(status, name, AT_);
1800
1801 // Extract strings from buffer
1802 for(size_type j = start_[0]; j < totalCount; j++) {
1803 MString tmp;
1804 tmp.append(&buf[(j - start_[0]) * strLen], strLen);
1805 data[j].append(tmp, 0, strLen);
1806 }
1807 }
1808 }
1809
1810 // Unpack strided data if necessary
1811 if(memoryStride != 1) {
1812 for(MPI_Offset i = 0; i < totalCount; i++) {
1813 array[memoryStride * i] = tmpScratch[i];
1814 }
1815 }
1816}
void b_ncEndDef()
Leave the define mode (NetCDF-specific). [MPI]

◆ b_readArray() [2/3]

template<class T >
void ParallelIoPNetcdf::b_readArray ( T *  array,
const MString path,
const MString name,
const size_type  noDims,
const size_type start,
const size_type count 
)
private

Definition at line 1595 of file parallelio_pnetcdf.cpp.

1596 {
1597 TRACE();
1598
1599 (void)array;
1600 (void)path;
1601 (void)name;
1602 (void)noDims;
1603 (void)start;
1604 (void)count;
1605 mTerm(1, AT_, "Group functionality not supported by PNetcdf backend");
1606}

◆ b_readArray() [3/3]

template<class T >
void ParallelIoPNetcdf::b_readArray ( T *  array,
const MString name,
const size_type  noDims,
MPI_Offset *  start,
MPI_Offset *  count,
MPI_Offset  memoryStride,
const size_type  noChunks,
MPI_Offset  diskStride 
)
private
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]arrayPointer to data in memory.
[in]nameName of the dataset.
[in]noDimsNumber of array dimensions
[in]startStart offset in the file (i.e. the domain should start reading from here).
[in]countNumber of elements read from this domain.
[in]memoryStrideStride in memory between data values.

This method is the exact counterpart to h5WriteArray().

Definition at line 1626 of file parallelio_pnetcdf.cpp.

1628 {
1629 TRACE();
1630
1631 b_ncEndDef();
1632
1633 // Get variable id
1634 MInt varId;
1635 MInt status = ncmpi_inq_varid(b_ncId, name.c_str(), &varId);
1636 b_error(status, name, AT_);
1637
1638 // Determine total data count
1639 size_type totalCount = 1;
1640 for(size_type d = 0; d < noDims; d++) {
1641 totalCount *= count[d];
1642 }
1643
1644 // Create temporary storage space if needed and set data pointers
1645 MInt tmpScratchSize = (memoryStride == 1 ? 1 : totalCount);
1646 ScratchSpace<T> tmpScratch(tmpScratchSize, FUN_, "tmpStorage");
1647 T* data = nullptr;
1648 if(memoryStride == 1) {
1649 data = array;
1650 } else {
1651 data = tmpScratch.data();
1652 }
1653
1654 // See b_writeArray(T*) for explanation
1655 if(count[0] == 0) {
1656 start[0] = 0;
1657 }
1658
1659 // Read array
1660 if(noChunks == 1) {
1661 // If number of chunks is one, read everything at once
1662 if(diskStride == 1) {
1663 status = pnetcdf_traits<T>::ncmpi_get_vara_type_all(b_ncId, varId, start, count, data);
1664 } else {
1665 status = pnetcdf_traits<T>::ncmpi_get_vars_type_all(b_ncId, varId, start, count, &diskStride, data);
1666 }
1667 b_error(status, name, AT_);
1668 } else {
1669 // Read in chunks
1670 MPI_Offset chunkSize = count[0] / noChunks;
1671 if(count[0] % noChunks > 0) {
1672 chunkSize += 1;
1673 }
1674
1675 // Determine number of entries for a fixed first dimension index
1676 size_type nDSize = 1;
1677 for(size_type d = 1; d < noDims; d++) {
1678 nDSize *= count[d];
1679 }
1680
1681 ScratchSpace<MPI_Offset> start_(noDims, FUN_, "start_");
1682 ScratchSpace<MPI_Offset> count_(noDims, FUN_, "count_");
1683
1684 std::copy(start, start + noDims, &start_[0]);
1685 std::copy(count, count + noDims, &count_[0]);
1686
1687 for(size_type i = 0; i < noChunks; i++) {
1688 start_[0] = min(start[0] + i * chunkSize * diskStride, count[0] * diskStride - 1);
1689 count_[0] = max(min(chunkSize, count[0] - i * chunkSize), 0ll);
1690 T* data_ = data + min(i * chunkSize * nDSize, (count[0] - 1) * nDSize);
1691 if(diskStride == 1) {
1692 status = pnetcdf_traits<T>::ncmpi_get_vara_type_all(b_ncId, varId, &start_[0], &count_[0], data_);
1693 } else {
1694 status = pnetcdf_traits<T>::ncmpi_get_vars_type_all(b_ncId, varId, &start_[0], &count_[0], &diskStride, data_);
1695 }
1696 b_error(status, name, AT_);
1697 }
1698 }
1699
1700 // Unpack strided data if necessary
1701 if(memoryStride != 1) {
1702 for(MPI_Offset i = 0; i < totalCount; i++) {
1703 array[memoryStride * i] = tmpScratch[i];
1704 }
1705 }
1706}

◆ b_readScalar() [1/2]

template<>
void ParallelIoPNetcdf::b_readScalar ( MString 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
Parameters
[out]scalarValue to read.
[in]nameName of the dataset.

This method is the exact counterpart to h5WriteScalar().

Definition at line 1864 of file parallelio_pnetcdf.cpp.

1864 {
1865 TRACE();
1866
1867 b_ncEndDef();
1868 MInt status;
1869
1870 // Get variable id
1871 MInt varId;
1872 status = ncmpi_inq_varid(b_ncId, name.c_str(), &varId);
1873 b_error(status, name, AT_);
1874
1875 // Determine offsets
1876 MPI_Offset start = 0;
1877 MPI_Offset count = 1;
1878
1879 // Read scalar (one char)
1880 status = ncmpi_get_vara_text_all(b_ncId, varId, &start, &count, (char*)scalar);
1881 b_error(status, name, AT_);
1882}

◆ b_readScalar() [2/2]

template<class T >
void ParallelIoPNetcdf::b_readScalar ( T *  scalar,
const MString name 
)
private
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
[out]scalarValue to read.
[in]nameName of the dataset.

This method is the exact counterpart to h5WriteScalar().

Definition at line 1831 of file parallelio_pnetcdf.cpp.

1831 {
1832 TRACE();
1833
1834 b_ncEndDef();
1835 MInt status;
1836
1837 // Get variable id
1838 MInt varId;
1839 status = ncmpi_inq_varid(b_ncId, name.c_str(), &varId);
1840 b_error(status, name, AT_);
1841
1842 // Determine offsets
1843 MPI_Offset start = 0;
1844 MPI_Offset count = 1;
1845
1846 // Read scalar
1847 status = pnetcdf_traits<T>::ncmpi_get_vara_type_all(b_ncId, varId, &start, &count, scalar);
1848 b_error(status, name, AT_);
1849}

◆ b_saveHeader()

void ParallelIoPNetcdf::b_saveHeader ( )
private
Author
Michael Schlottke (mic) mic@a.nosp@m.ia.r.nosp@m.wth-a.nosp@m.ache.nosp@m.n.de
Date
2013-04-09

This method must be called exactly once for newly created files.

Definition at line 538 of file parallelio_pnetcdf.cpp.

538 {
539 TRACE();
540
543}
void b_writeAdditionalData()
Write additional data to file (NetCDF-specific). [MPI]
void b_addAdditionalHeader()
Write additional headers to file (e.g. grid file name, creation date etc.). [MPI]

◆ b_setAttribute() [1/2]

template<>
void ParallelIoPNetcdf::b_setAttribute ( const MString value,
const MString name,
const MString datasetName,
const size_type  totalCount 
)
Author
Ramandeep Jain (HiWi) raman.nosp@m.deep.nosp@m.jain@.nosp@m.gmai.nosp@m.l.com, Konstantin Froehlich
Date
2015-Mar
Parameters
[in]valueAttribute value.
[in]nameAttribute name.
[in]datasetNameDataset (scalar or array) to which the attribute is attached (leave empty for file attributes).

Definition at line 1938 of file parallelio_pnetcdf.cpp.

1939 {
1940 TRACE();
1941
1942 if(totalCount > 1) {
1943 mTerm(1, AT_, "Array of strings attributes not yet supported.");
1944 }
1945
1946 // Determine variable id
1947 MInt varId;
1948 if(datasetName.empty()) {
1949 varId = NC_GLOBAL;
1950 } else {
1951 // If this is a data file, get data file-specific name
1952 MInt status = ncmpi_inq_varid(b_ncId, datasetName.c_str(), &varId);
1953 b_error(status, datasetName, AT_);
1954 }
1955
1956 // If attribute does not exist or is of greater size, go to define mode
1957 if(!b_hasAttribute(name, datasetName)) {
1958 b_ncRedef();
1959 } else {
1960 MPI_Offset length;
1961 MInt status = ncmpi_inq_attlen(b_ncId, varId, name.c_str(), &length);
1962 b_error(status, datasetName + "::" + name, AT_);
1963
1964 if(length < static_cast<MPI_Offset>(value->length())) {
1965 b_ncRedef();
1966 }
1967 }
1968
1969 // Write attribute
1970 MInt status = ncmpi_put_att_text(b_ncId, varId, name.c_str(), value->length(), (*value).c_str());
1971 b_error(status, datasetName + "::" + name, AT_);
1972}
MBool b_hasAttribute(const MString &name, const MString &datasetName="")
Check if a given attribute exists in the file.

◆ b_setAttribute() [2/2]

template<class T >
void ParallelIoPNetcdf::b_setAttribute ( const T *  value,
const MString name,
const MString datasetName = "",
const size_type  totalCount = 1 
)
private
Author
Ramandeep Jain (HiWi) raman.nosp@m.deep.nosp@m.jain@.nosp@m.gmai.nosp@m.l.com, Konstantin Froehlich
Date
2015-Mar1
Parameters
[in]valueAttribute value.
[in]nameAttribute name.
[in]datasetNameDataset (scalar or array) to which the attribute is attached (leave empty for file attributes).

Definition at line 1901 of file parallelio_pnetcdf.cpp.

1902 {
1903 TRACE();
1904
1905 // Determine variable id
1906 MInt varId;
1907 if(datasetName.empty()) {
1908 varId = NC_GLOBAL;
1909 } else {
1910 MInt status = ncmpi_inq_varid(b_ncId, datasetName.c_str(), &varId);
1911 b_error(status, datasetName, AT_);
1912 }
1913
1914 // If attribute does not exist, go to define mode
1915 if(!b_hasAttribute(name, datasetName)) {
1916 b_ncRedef();
1917 }
1918
1919 // Write attribute
1920 MInt status =
1921 pnetcdf_traits<T>::ncmpi_put_att_type(b_ncId, varId, name.c_str(), pnetcdf_traits<T>::type(), totalCount, value);
1922 b_error(status, datasetName + "::" + name, AT_);
1923}

◆ b_writeAdditionalData()

void ParallelIoPNetcdf::b_writeAdditionalData ( )
private
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

Definition at line 683 of file parallelio_pnetcdf.cpp.

683 {
684 TRACE();
685
686 // At the moment, nothing happens here
687}

◆ b_writeArray() [1/3]

template<>
void ParallelIoPNetcdf::b_writeArray ( const MString array,
const MString name,
const size_type  noDims,
MPI_Offset *  start,
MPI_Offset *  count,
MPI_Offset  memoryStride,
const size_type  noChunks,
MPI_Offset  diskStride 
)
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]arrayPointer to data in memory.
[in]nameName of the dataset.
[in]noDimsNumber of array dimensions
[in]startStart offset in the file (i.e. the domain should start writing from here).
[in]countNumber of elements written from this domain.
[in]memoryStrideStride in memory between data values.

Definition at line 1436 of file parallelio_pnetcdf.cpp.

1438 {
1439 TRACE();
1440
1441 b_ncEndDef();
1442
1443 // Get variable id
1444 MInt varId;
1445 MInt status = ncmpi_inq_varid(b_ncId, name.c_str(), &varId);
1446 b_error(status, name, AT_);
1447
1448 // Determine total number of strings (last dimension is string length)
1449 size_type totalCount = 1;
1450 for(size_type d = 0; d < noDims - 1; d++) {
1451 totalCount *= count[d];
1452 }
1453
1454 // Determine length of one string
1455 size_type strLen = count[noDims - 1];
1456
1457 // Create temporary storage space if needed and set data pointers
1458 MInt tmpScratchSize = (memoryStride == 1 ? 1 : totalCount);
1459 ScratchSpace<MString> tmpScratch(tmpScratchSize, FUN_, "tmpStorage");
1460
1461 // Pack strided data
1462 const MString* data = 0;
1463 if(memoryStride == 1) {
1464 data = array;
1465 } else {
1466 for(MPI_Offset i = 0; i < totalCount; i++) {
1467 tmpScratch[i] = array[memoryStride * i];
1468 }
1469 data = tmpScratch.data();
1470 }
1471
1472 // labels:IO this is a bugfix for the case when the last process in the communicator
1473 // has zero elements to write, which resulted in the error:
1474 //"NetCDF returns status -40: Index exceeds dimension bound"
1475 // solution taken from here:
1476 // http://lists.mcs.anl.gov/pipermail/parallel-netcdf/2004-December/000388.html
1477 if(count[0] == 0) {
1478 start[0] = 0;
1479 }
1480
1481 // Create buffer for writing
1482 size_type nCount = totalCount * strLen;
1483 ScratchSpace<char> buf(nCount, FUN_, "buf");
1484
1485 // Copy data to buffer
1486 for(MPI_Offset i = 0; i < totalCount; i++) {
1487 data[i].copy(&buf[i * strLen], strLen, 0);
1488 }
1489
1490 // Write array
1491 if(noChunks == 1) {
1492 // If number of chunks is one, write everything at once
1493 if(diskStride == 1) {
1494 status = ncmpi_put_vara_text_all(b_ncId, varId, start, count, &buf[0]);
1495 } else {
1496 status = ncmpi_put_vars_text_all(b_ncId, varId, start, count, &diskStride, &buf[0]);
1497 }
1498 b_error(status, name, AT_);
1499 } else {
1500 // Write in chunks
1501 MPI_Offset chunkSize = count[0] / noChunks;
1502 if(count[0] % noChunks > 0) {
1503 chunkSize += 1;
1504 }
1505
1506 // Determine number of entries for a fixed first dimension index
1507 size_type nDSize = 1;
1508 for(size_type d = 1; d < noDims; d++) {
1509 nDSize *= count[d];
1510 }
1511
1512 ScratchSpace<MPI_Offset> start_(noDims, FUN_, "start_");
1513 ScratchSpace<MPI_Offset> count_(noDims, FUN_, "count_");
1514
1515 std::copy(start, start + noDims, &start_[0]);
1516 std::copy(count, count + noDims, &count_[0]);
1517
1518 for(size_type i = 0; i < noChunks; i++) {
1519 start_[0] = min(start[0] + i * chunkSize * diskStride, count[0] * diskStride - 1);
1520 count_[0] = max(min(chunkSize, count[0] - i * chunkSize), 0ll);
1521 const char* buf_ = &buf[0] + min(i * chunkSize * nDSize, (count[0] - 1) * nDSize);
1522 if(diskStride == 1) {
1523 status = ncmpi_put_vara_text_all(b_ncId, varId, &start_[0], &count_[0], &buf_[0]);
1524 } else {
1525 status = ncmpi_put_vars_text_all(b_ncId, varId, &start_[0], &count_[0], &diskStride, &buf_[0]);
1526 }
1527 b_error(status, name, AT_);
1528 }
1529 }
1530}

◆ b_writeArray() [2/3]

template<class T >
void ParallelIoPNetcdf::b_writeArray ( const T *  array,
const MString path,
const MString name,
const size_type  noDims,
const size_type start,
const size_type count,
const size_type ghost 
)
private

Definition at line 1303 of file parallelio_pnetcdf.cpp.

1304 {
1305 (void)array;
1306 (void)path;
1307 (void)name;
1308 (void)noDims;
1309 (void)start;
1310 (void)count;
1311 (void)ghost;
1312 TRACE();
1313 mTerm(1, AT_, "Group functionality not supported by PNetcdf backend");
1314}

◆ b_writeArray() [3/3]

template<class T >
void ParallelIoPNetcdf::b_writeArray ( const T *const  array,
const MString name,
const size_type  noDims,
MPI_Offset *  start,
MPI_Offset *  count,
MPI_Offset  memoryStride,
const size_type  noChunks,
MPI_Offset  diskStride 
)
private
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 data to write.
Parameters
[in]arrayPointer to data in memory.
[in]nameName of the dataset.
[in]noDimsNumber of array dimensions
[in]startStart offset in the file (i.e. the domain should start writing from here).
[in]countNumber of elements written from this domain.
[in]memoryStrideStride in memory between data values.

Definition at line 1336 of file parallelio_pnetcdf.cpp.

1338 {
1339 TRACE();
1340
1341 b_ncEndDef();
1342
1343 // Get variable id
1344 MInt varId;
1345 MInt status = ncmpi_inq_varid(b_ncId, name.c_str(), &varId);
1346 b_error(status, name, AT_);
1347
1348 // Determine total data count
1349 size_type totalCount = 1;
1350 for(size_type d = 0; d < noDims; d++) {
1351 totalCount *= count[d];
1352 }
1353
1354 // Create temporary storage space if needed and set data pointers
1355 MInt tmpScratchSize = (memoryStride == 1) ? 1 : totalCount;
1356 ScratchSpace<T> tmpScratch(tmpScratchSize, FUN_, "tmpStorage");
1357
1358 // Pack strided data
1359 const T* data = 0;
1360 if(memoryStride == 1) {
1361 data = array;
1362 } else {
1363 for(MPI_Offset i = 0; i < totalCount; i++) {
1364 tmpScratch[i] = array[memoryStride * i];
1365 }
1366 data = tmpScratch.data();
1367 }
1368
1369 // labels:IO this is a bugfix for the case when the last process in the communicator
1370 // has zero elements to write, which resulted in the error:
1371 //"NetCDF returns status -40: Index exceeds dimension bound"
1372 // solution taken from here:
1373 // http://lists.mcs.anl.gov/pipermail/parallel-netcdf/2004-December/000388.html
1374 if(count[0] == 0) {
1375 start[0] = 0;
1376 }
1377
1378 // Write array
1379 if(noChunks == 1) {
1380 // If number of chunks is one, write everything at once
1381 if(diskStride == 1) {
1382 status = pnetcdf_traits<T>::ncmpi_put_vara_type_all(b_ncId, varId, start, count, data);
1383 } else {
1384 status = pnetcdf_traits<T>::ncmpi_put_vars_type_all(b_ncId, varId, start, count, &diskStride, data);
1385 }
1386 b_error(status, name, AT_);
1387 } else {
1388 // Write in chunks
1389 MPI_Offset chunkSize = count[0] / noChunks;
1390 if(count[0] % noChunks > 0) {
1391 chunkSize += 1;
1392 }
1393
1394 // Determine number of entries for a fixed first dimension index
1395 size_type nDSize = 1;
1396 for(size_type d = 1; d < noDims; d++) {
1397 nDSize *= count[d];
1398 }
1399
1400 ScratchSpace<MPI_Offset> start_(noDims, FUN_, "start_");
1401 ScratchSpace<MPI_Offset> count_(noDims, FUN_, "count_");
1402
1403 std::copy(start, start + noDims, &start_[0]);
1404 std::copy(count, count + noDims, &count_[0]);
1405
1406 for(size_type i = 0; i < noChunks; i++) {
1407 start_[0] = min(start[0] + i * chunkSize * diskStride, count[0] * diskStride - 1);
1408 count_[0] = max(min(chunkSize, count[0] - i * chunkSize), 0ll);
1409 const T* data_ = data + min(i * chunkSize * nDSize, (count[0] - 1) * nDSize);
1410 if(diskStride == 1) {
1411 status = pnetcdf_traits<T>::ncmpi_put_vara_type_all(b_ncId, varId, &start_[0], &count_[0], data_);
1412 } else {
1413 status = pnetcdf_traits<T>::ncmpi_put_vars_type_all(b_ncId, varId, &start_[0], &count_[0], &diskStride, data_);
1414 }
1415 b_error(status, name, AT_);
1416 }
1417 }
1418}

◆ b_writeScalar() [1/2]

template<>
void ParallelIoPNetcdf::b_writeScalar ( const MString 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
Parameters
[in]scalarValue to write.
[in]nameName of the dataset.

Definition at line 1574 of file parallelio_pnetcdf.cpp.

1574 {
1575 TRACE();
1576
1577 b_ncEndDef();
1578 MInt status;
1579
1580 // Get variable id
1581 MInt varId;
1582 status = ncmpi_inq_varid(b_ncId, name.c_str(), &varId);
1583 b_error(status, name, AT_);
1584
1585 // Determine offsets
1586 MPI_Offset start = 0;
1587 MPI_Offset count = 1;
1588
1589 // Write scalar
1590 status = ncmpi_put_vara_text_all(b_ncId, varId, &start, &count, scalar.c_str());
1591 b_error(status, name, AT_);
1592}

◆ b_writeScalar() [2/2]

template<class T >
void ParallelIoPNetcdf::b_writeScalar ( const T  scalar,
const MString name 
)
private
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]scalarValue to write.
[in]nameName of the dataset.

Definition at line 1543 of file parallelio_pnetcdf.cpp.

1543 {
1544 TRACE();
1545
1546 b_ncEndDef();
1547 MInt status;
1548
1549 // Get variable id
1550 MInt varId;
1551 status = ncmpi_inq_varid(b_ncId, name.c_str(), &varId);
1552 b_error(status, name, AT_);
1553
1554 // Determine offsets
1555 MPI_Offset start = 0;
1556 MPI_Offset count = 1;
1557
1558 // Write scalar
1559 status = pnetcdf_traits<T>::ncmpi_put_vara_type_all(b_ncId, varId, &start, &count, &scalar);
1560 b_error(status, name, AT_);
1561}

◆ close()

void ParallelIoPNetcdf::close ( )

Definition at line 465 of file parallelio_pnetcdf.cpp.

465 {
466 TRACE();
467#ifdef MAIA_NCMPI_PRINT_FILE_HINTS
468 if(m_domainId == 0) {
469 std::cerr << std::endl << "Closing file: " << m_fileName << std::endl;
471 }
472#endif
473
474 MInt status = ncmpi_close(b_ncId);
475 b_error(status, m_fileName, AT_);
476 b_ncId = -1; // Reset Netcdf file id, calling any other function after close() should result in an invalid file id
477 // Netcdf error
478
479#ifdef MAIA_EXTRA_DEBUG
480 for(set<MString>::iterator it = m_unwrittenArrays.begin(); it != m_unwrittenArrays.end(); ++it) {
481 cerr << "Warning: array '" << *it << "' in file '" << m_fileName << "' was defined but never written. "
482 << "Make sure that this is the intended behavior." << endl;
483 }
484 for(set<MString>::iterator it = m_unwrittenScalars.begin(); it != m_unwrittenScalars.end(); ++it) {
485 cerr << "Warning: scalar '" << *it << "' in file '" << m_fileName << "' was defined but never written. "
486 << "Make sure that this is the intended behavior." << endl;
487 }
488#endif
489}
std::set< MString > m_unwrittenScalars
Definition: parallelio.h:255
std::set< MString > m_unwrittenArrays
Definition: parallelio.h:254

Member Data Documentation

◆ b_ncDataMode

MBool ParallelIoPNetcdf::b_ncDataMode = false
private

Definition at line 92 of file parallelio_pnetcdf.h.

◆ b_ncDimensions

NcDimMap ParallelIoPNetcdf::b_ncDimensions {}
private

Definition at line 97 of file parallelio_pnetcdf.h.

◆ b_ncId

MInt ParallelIoPNetcdf::b_ncId = -1
private

Definition at line 91 of file parallelio_pnetcdf.h.

◆ ParallelIoBase< ParallelIoPNetcdf >

Definition at line 13 of file parallelio_pnetcdf.h.


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