Nektar++
Classes | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Private Member Functions | List of all members
Nektar::LibUtilities::FieldIOHdf5 Class Reference

#include <FieldIOHdf5.h>

Inheritance diagram for Nektar::LibUtilities::FieldIOHdf5:
[legend]

Classes

struct  OffsetHelper
 

Public Member Functions

 FieldIOHdf5 (LibUtilities::CommSharedPtr pComm, bool sharedFilesystem)
 Construct the FieldIO object for HDF5 output. More...
 
 ~FieldIOHdf5 () override
 
- Public Member Functions inherited from Nektar::LibUtilities::FieldIO
 FieldIO (LibUtilities::CommSharedPtr pComm, bool sharedFilesystem)
 Constructor for FieldIO base class. More...
 
virtual ~FieldIO ()
 
void Write (const std::string &outFile, std::vector< FieldDefinitionsSharedPtr > &fielddefs, std::vector< std::vector< NekDouble > > &fielddata, const FieldMetaDataMap &fieldinfomap=NullFieldMetaDataMap, const bool backup=false)
 Write out the field information to the file outFile. More...
 
void Import (const std::string &infilename, std::vector< FieldDefinitionsSharedPtr > &fielddefs, std::vector< std::vector< NekDouble > > &fielddata=NullVectorNekDoubleVector, FieldMetaDataMap &fieldinfomap=NullFieldMetaDataMap, const Array< OneD, int > &ElementIDs=NullInt1DArray)
 Read field information from the file infilename. More...
 
DataSourceSharedPtr ImportFieldMetaData (const std::string &filename, FieldMetaDataMap &fieldmetadatamap)
 Import the metadata from a field file. More...
 
const std::string & GetClassName () const
 
std::string GetFileEnding () const
 Helper function that determines default file extension. More...
 

Static Public Member Functions

static FieldIOSharedPtr create (LibUtilities::CommSharedPtr pComm, bool sharedFilesystem)
 Creates an instance of this class. More...
 
- Static Public Member Functions inherited from Nektar::LibUtilities::FieldIO
static const std::string GetFileType (const std::string &filename, CommSharedPtr comm)
 Determine file type of given input file. More...
 
static std::shared_ptr< FieldIOCreateDefault (const LibUtilities::SessionReaderSharedPtr session)
 Returns an object for the default FieldIO method. More...
 
static std::shared_ptr< FieldIOCreateForFile (const LibUtilities::SessionReaderSharedPtr session, const std::string &filename)
 Construct a FieldIO object for a given input filename. More...
 
static void AddInfoTag (TagWriterSharedPtr root, const FieldMetaDataMap &fieldmetadatamap)
 Add provenance information to the field metadata map. More...
 

Static Public Attributes

static const unsigned int FORMAT_VERSION = 1
 Version of the Nektar++ HDF5 format, which is embedded into the main NEKTAR group as an attribute. More...
 
static const unsigned int ELEM_DCMP_IDX = 0
 A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of elements in decomposition (i.e. field definition). More...
 
static const unsigned int VAL_DCMP_IDX = 1
 A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of data points in decomposition (i.e. field definition). More...
 
static const unsigned int ORDER_DCMP_IDX = 2
 A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of elements multiplied by the dimension of the element, giving number of modes when variable polynomial order is defined. More...
 
static const unsigned int HOMY_DCMP_IDX = 3
 A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of the number of y-planes for homogeneous simulations. More...
 
static const unsigned int HOMZ_DCMP_IDX = 4
 A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of the number of z-planes for homogeneous simulations. More...
 
static const unsigned int HOMS_DCMP_IDX = 5
 A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of the number of strips for homogeneous simulations. More...
 
static const unsigned int HASH_DCMP_IDX = 6
 The hash of the field definition information, which defines the name of the attribute containing the field definition itself. More...
 
static const unsigned int MAX_DCMPS = FieldIOHdf5::HASH_DCMP_IDX + 1
 A helper for FieldIOHdf5::v_Write. Describes the maximum number of items in the decomposition per field definition. More...
 
static const unsigned int ELEM_CNT_IDX = 0
 A helper for FieldIOHdf5::v_Write. Describes the position of the number of elements in the cnt array. More...
 
static const unsigned int VAL_CNT_IDX = 1
 A helper for FieldIOHdf5::v_Write. Describes the position of the number of data points in the cnt array. More...
 
static const unsigned int ORDER_CNT_IDX = 2
 A helper for FieldIOHdf5::v_Write. Describes the position of the number of order points in the cnt array. More...
 
static const unsigned int HOMY_CNT_IDX = 3
 A helper for FieldIOHdf5::v_Write. Describes the position of the number of homogeneous y-planes in the cnt array. More...
 
static const unsigned int HOMZ_CNT_IDX = 4
 A helper for FieldIOHdf5::v_Write. Describes the position of the number of homogeneous z-planes in the cnt array. More...
 
static const unsigned int HOMS_CNT_IDX = 5
 A helper for FieldIOHdf5::v_Write. Describes the position of the number of homogeneous strips in the cnt array. More...
 
static const unsigned int MAX_CNTS = FieldIOHdf5::HOMS_CNT_IDX + 1
 A helper for FieldIOHdf5::v_Write. Describes the maximum number of items in the cnt array per field definition. More...
 
static const unsigned int IDS_IDX_IDX = 0
 A helper for FieldIOHdf5::v_Write. Describes the position of the element IDs within the indexing set. More...
 
static const unsigned int DATA_IDX_IDX = 1
 A helper for FieldIOHdf5::v_Write. Describes the position of the data size within the indexing set. More...
 
static const unsigned int ORDER_IDX_IDX = 2
 A helper for FieldIOHdf5::v_Write. Describes the position of the element order within the indexing set. More...
 
static const unsigned int HOMY_IDX_IDX = 3
 A helper for FieldIOHdf5::v_Write. Describes the position of the number of y-planes within the indexing set. More...
 
static const unsigned int HOMZ_IDX_IDX = 4
 A helper for FieldIOHdf5::v_Write. Describes the position of the number of z-planes within the indexing set. More...
 
static const unsigned int HOMS_IDX_IDX = 5
 A helper for FieldIOHdf5::v_Write. Describes the position of the number of homogeneous strips within the indexing set. More...
 
static const unsigned int MAX_IDXS = FieldIOHdf5::HOMS_IDX_IDX + 1
 A helper for FieldIOHdf5::v_Write. Describes the maximum number of items in the indexing set. More...
 
static std::string className
 Name of class. More...
 

Protected Member Functions

void v_Write (const std::string &outFile, std::vector< FieldDefinitionsSharedPtr > &fielddefs, std::vector< std::vector< NekDouble > > &fielddata, const FieldMetaDataMap &fieldinfomap=NullFieldMetaDataMap, const bool backup=false) override
 Write a HDF5 file to outFile given the field definitions fielddefs, field data fielddata and metadata fieldmetadatamap. More...
 
void v_Import (const std::string &infilename, std::vector< FieldDefinitionsSharedPtr > &fielddefs, std::vector< std::vector< NekDouble > > &fielddata=NullVectorNekDoubleVector, FieldMetaDataMap &fieldinfomap=NullFieldMetaDataMap, const Array< OneD, int > &ElementIDs=NullInt1DArray) override
 Import a HDF5 format file. More...
 
DataSourceSharedPtr v_ImportFieldMetaData (const std::string &filename, FieldMetaDataMap &fieldmetadatamap) override
 Import field metadata from filename and return the data source which wraps filename. More...
 
const std::string & v_GetClassName () const override
 Get class name. More...
 
- Protected Member Functions inherited from Nektar::LibUtilities::FieldIO
int CheckFieldDefinition (const FieldDefinitionsSharedPtr &fielddefs)
 Check field definitions for correctness and return storage size. More...
 
std::vector< unsigned int > GetNumberOfCoeffsPerElement (const FieldDefinitionsSharedPtr &fielddefs)
 Compute number of data points needed to store expansion inside each element. More...
 
std::string SetUpOutput (const std::string outname, bool perRank, bool backup=false)
 Set up the filesystem ready for output. More...
 
virtual void v_Write (const std::string &outFile, std::vector< FieldDefinitionsSharedPtr > &fielddefs, std::vector< std::vector< NekDouble > > &fielddata, const FieldMetaDataMap &fieldinfomap, const bool backup=false)=0
 Write out the field information to the file outFile. More...
 
virtual void v_Import (const std::string &infilename, std::vector< FieldDefinitionsSharedPtr > &fielddefs, std::vector< std::vector< NekDouble > > &fielddata=NullVectorNekDoubleVector, FieldMetaDataMap &fieldinfomap=NullFieldMetaDataMap, const Array< OneD, int > &ElementIDs=NullInt1DArray)=0
 Read field information from the file infilename. More...
 
virtual DataSourceSharedPtr v_ImportFieldMetaData (const std::string &filename, FieldMetaDataMap &fieldmetadatamap)=0
 Import the metadata from a field file. More...
 
virtual std::string v_GetFileEnding () const
 
virtual const std::string & v_GetClassName () const =0
 

Private Member Functions

void ImportHDF5FieldMetaData (DataSourceSharedPtr dataSource, FieldMetaDataMap &fieldmetadatamap)
 Import field metadata from dataSource. More...
 
void ImportFieldDef (H5::PListSharedPtr readPL, H5::GroupSharedPtr root, std::vector< uint64_t > &decomps, uint64_t decomp, OffsetHelper offset, std::string group, FieldDefinitionsSharedPtr def)
 Import field definitions from a HDF5 document. More...
 

Additional Inherited Members

- Protected Attributes inherited from Nektar::LibUtilities::FieldIO
LibUtilities::CommSharedPtr m_comm
 Communicator to use when writing parallel format. More...
 
bool m_sharedFilesystem
 Boolean dictating whether we are on a shared filesystem. More...
 

Detailed Description

Definition at line 173 of file FieldIOHdf5.h.

Constructor & Destructor Documentation

◆ FieldIOHdf5()

Nektar::LibUtilities::FieldIOHdf5::FieldIOHdf5 ( LibUtilities::CommSharedPtr  pComm,
bool  sharedFilesystem 
)

Construct the FieldIO object for HDF5 output.

Parameters
pCommCommunicator object.
sharedFilesystemTrue if this system has a shared filesystem.

Definition at line 143 of file FieldIOHdf5.cpp.

145 : FieldIO(pComm, sharedFilesystem)
146{
147}
FieldIO(LibUtilities::CommSharedPtr pComm, bool sharedFilesystem)
Constructor for FieldIO base class.
Definition: FieldIO.cpp:321

◆ ~FieldIOHdf5()

Nektar::LibUtilities::FieldIOHdf5::~FieldIOHdf5 ( )
inlineoverride

Definition at line 217 of file FieldIOHdf5.h.

218 {
219 }

Member Function Documentation

◆ create()

static FieldIOSharedPtr Nektar::LibUtilities::FieldIOHdf5::create ( LibUtilities::CommSharedPtr  pComm,
bool  sharedFilesystem 
)
inlinestatic

Creates an instance of this class.

Definition at line 204 of file FieldIOHdf5.h.

206 {
208 sharedFilesystem);
209 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

◆ ImportFieldDef()

void Nektar::LibUtilities::FieldIOHdf5::ImportFieldDef ( H5::PListSharedPtr  readPL,
H5::GroupSharedPtr  root,
std::vector< uint64_t > &  decomps,
uint64_t  decomp,
OffsetHelper  offset,
std::string  group,
FieldDefinitionsSharedPtr  def 
)
private

Import field definitions from a HDF5 document.

Parameters
readPLReading parameter list.
rootRoot group containing field definitions.
groupGroup name to process.
defOn output contains field definitions.

Definition at line 1191 of file FieldIOHdf5.cpp.

1197{
1198 std::stringstream prfx;
1199 prfx << m_comm->GetRank() << ": FieldIOHdf5::ImportFieldDefsHdf5(): ";
1200
1201 H5::GroupSharedPtr field = root->OpenGroup(group);
1202 ASSERTL1(field, prfx.str() + "cannot open field group, " + group + '.');
1203
1204 def->m_uniOrder = false;
1205
1206 H5::Group::AttrIterator attrIt = field->attr_begin();
1207 H5::Group::AttrIterator attrEnd = field->attr_end();
1208 for (; attrIt != attrEnd; ++attrIt)
1209 {
1210 const std::string &attrName = *attrIt;
1211 if (attrName == "FIELDS")
1212 {
1213 field->GetAttribute(attrName, def->m_fields);
1214 }
1215 else if (attrName == "SHAPE")
1216 {
1217 std::string shapeString;
1218 field->GetAttribute(attrName, shapeString);
1219
1220 // check to see if homogeneous expansion and if so
1221 // strip down the shapeString definition
1222 size_t loc;
1223 //---> this finds the first location of 'n'!
1224 if (shapeString.find("Strips") != std::string::npos)
1225 {
1226 def->m_homoStrips = true;
1227 }
1228
1229 if ((loc = shapeString.find_first_of("-")) != std::string::npos)
1230 {
1231 if (shapeString.find("Exp1D") != std::string::npos)
1232 {
1233 def->m_numHomogeneousDir = 1;
1234 }
1235 else // HomogeneousExp1D
1236 {
1237 def->m_numHomogeneousDir = 2;
1238 }
1239
1240 shapeString.erase(loc, shapeString.length());
1241 }
1242
1243 // get the geometrical shape
1244 bool valid = false;
1245 for (unsigned int j = 0; j < SIZE_ShapeType; j++)
1246 {
1247 if (ShapeTypeMap[j] == shapeString)
1248 {
1249 def->m_shapeType = (ShapeType)j;
1250 valid = true;
1251 break;
1252 }
1253 }
1254
1255 ASSERTL0(
1256 valid,
1257 prfx.str() +
1258 std::string("unable to correctly parse the shape type: ")
1259 .append(shapeString)
1260 .c_str());
1261 }
1262 else if (attrName == "BASIS")
1263 {
1264 field->GetAttribute(attrName, def->m_basis);
1265 // check the basis is in range
1266 for (auto &bIt : def->m_basis)
1267 {
1268 BasisType bt = bIt;
1269 ASSERTL0(bt >= 0 && bt < SIZE_BasisType,
1270 prfx.str() +
1271 "unable to correctly parse the basis types.");
1272 }
1273 }
1274 else if (attrName == "HOMOGENEOUSLENGTHS")
1275 {
1276 field->GetAttribute(attrName, def->m_homogeneousLengths);
1277 }
1278 else if (attrName == "NUMMODESPERDIR")
1279 {
1280 std::string numModesPerDir;
1281 field->GetAttribute(attrName, numModesPerDir);
1282
1283 if (strstr(numModesPerDir.c_str(), "UNIORDER:"))
1284 {
1285 def->m_uniOrder = true;
1286 bool valid = ParseUtils::GenerateVector(
1287 numModesPerDir.substr(9), def->m_numModes);
1288 ASSERTL0(valid,
1289 prfx.str() +
1290 "unable to correctly parse the number of modes.");
1291 }
1292 }
1293 else if (attrName == "POINTSTYPE")
1294 {
1295 std::string pointsString;
1296 field->GetAttribute(attrName, pointsString);
1297 def->m_pointsDef = true;
1298
1299 std::vector<std::string> pointsStrings;
1300 bool valid =
1301 ParseUtils::GenerateVector(pointsString, pointsStrings);
1302 ASSERTL0(valid, prfx.str() +
1303 "unable to correctly parse the points types.");
1304 for (size_t i = 0; i < pointsStrings.size(); i++)
1305 {
1306 valid = false;
1307 for (unsigned int j = 0; j < SIZE_PointsType; j++)
1308 {
1309 if (kPointsTypeStr[j] == pointsStrings[i])
1310 {
1311 def->m_points.push_back((PointsType)j);
1312 valid = true;
1313 break;
1314 }
1315 }
1316
1317 ASSERTL0(valid,
1318 prfx.str() +
1319 std::string(
1320 "unable to correctly parse the points type: ")
1321 .append(pointsStrings[i])
1322 .c_str());
1323 }
1324 }
1325 else if (attrName == "NUMPOINTSPERDIR")
1326 {
1327 std::string numPointsPerDir;
1328 field->GetAttribute(attrName, numPointsPerDir);
1329 def->m_numPointsDef = true;
1330
1331 bool valid =
1332 ParseUtils::GenerateVector(numPointsPerDir, def->m_numPoints);
1333 ASSERTL0(valid,
1334 prfx.str() +
1335 "unable to correctly parse the number of points.");
1336 }
1337 else
1338 {
1339 std::string errstr("unknown attribute: ");
1340 errstr += attrName;
1341 ASSERTL1(false, prfx.str() + errstr.c_str());
1342 }
1343 }
1344
1345 if (def->m_numHomogeneousDir >= 1)
1346 {
1347 H5::DataSetSharedPtr homz_dset = root->OpenDataSet("HOMOGENEOUSZIDS");
1348 H5::DataSpaceSharedPtr homz_fspace = homz_dset->GetSpace();
1349 uint64_t dsize = decomps[decomp * MAX_DCMPS + HOMZ_DCMP_IDX];
1350 homz_fspace->SelectRange(offset.homz, dsize);
1351 homz_dset->Read(def->m_homogeneousZIDs, homz_fspace, readPL);
1352 }
1353
1354 if (def->m_numHomogeneousDir >= 2)
1355 {
1356 H5::DataSetSharedPtr homy_dset = root->OpenDataSet("HOMOGENEOUSYIDS");
1357 H5::DataSpaceSharedPtr homy_fspace = homy_dset->GetSpace();
1358 uint64_t dsize = decomps[decomp * MAX_DCMPS + HOMY_DCMP_IDX];
1359 homy_fspace->SelectRange(offset.homy, dsize);
1360 homy_dset->Read(def->m_homogeneousYIDs, homy_fspace, readPL);
1361 }
1362
1363 if (def->m_homoStrips)
1364 {
1365 H5::DataSetSharedPtr homs_dset = root->OpenDataSet("HOMOGENEOUSSIDS");
1366 H5::DataSpaceSharedPtr homs_fspace = homs_dset->GetSpace();
1367 uint64_t dsize = decomps[decomp * MAX_DCMPS + HOMS_DCMP_IDX];
1368 homs_fspace->SelectRange(offset.homs, dsize);
1369 homs_dset->Read(def->m_homogeneousSIDs, homs_fspace, readPL);
1370 }
1371
1372 if (!def->m_uniOrder)
1373 {
1374 H5::DataSetSharedPtr order_dset = root->OpenDataSet("POLYORDERS");
1375 H5::DataSpaceSharedPtr order_fspace = order_dset->GetSpace();
1376 uint64_t dsize = decomps[decomp * MAX_DCMPS + ORDER_DCMP_IDX];
1377 order_fspace->SelectRange(offset.order, dsize);
1378 order_dset->Read(def->m_numModes, order_fspace, readPL);
1379 }
1380}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
static const unsigned int ORDER_DCMP_IDX
A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of ...
Definition: FieldIOHdf5.h:180
static const unsigned int HOMS_DCMP_IDX
A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of ...
Definition: FieldIOHdf5.h:183
static const unsigned int HOMZ_DCMP_IDX
A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of ...
Definition: FieldIOHdf5.h:182
static const unsigned int MAX_DCMPS
A helper for FieldIOHdf5::v_Write. Describes the maximum number of items in the decomposition per fie...
Definition: FieldIOHdf5.h:185
static const unsigned int HOMY_DCMP_IDX
A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of ...
Definition: FieldIOHdf5.h:181
LibUtilities::CommSharedPtr m_comm
Communicator to use when writing parallel format.
Definition: FieldIO.h:280
static bool GenerateVector(const std::string &str, std::vector< T > &out)
Takes a comma-separated string and converts it to entries in a vector.
Definition: ParseUtils.cpp:130
std::shared_ptr< DataSpace > DataSpaceSharedPtr
Definition: H5.h:78
std::shared_ptr< Group > GroupSharedPtr
Definition: FieldIOHdf5.h:48
std::shared_ptr< DataSet > DataSetSharedPtr
Definition: H5.h:90
const char *const ShapeTypeMap[SIZE_ShapeType]
Definition: ShapeType.hpp:75
const std::string kPointsTypeStr[]
Definition: Foundations.hpp:52
@ SIZE_PointsType
Length of enum list.
Definition: PointsType.h:95
@ SIZE_BasisType
Length of enum list.
Definition: BasisType.h:70

References ASSERTL0, ASSERTL1, FilterPython_Function::field, Nektar::ParseUtils::GenerateVector(), Nektar::LibUtilities::FieldIOHdf5::OffsetHelper::homs, HOMS_DCMP_IDX, Nektar::LibUtilities::FieldIOHdf5::OffsetHelper::homy, HOMY_DCMP_IDX, Nektar::LibUtilities::FieldIOHdf5::OffsetHelper::homz, HOMZ_DCMP_IDX, Nektar::LibUtilities::kPointsTypeStr, CG_Iterations::loc, Nektar::LibUtilities::FieldIO::m_comm, MAX_DCMPS, Nektar::LibUtilities::FieldIOHdf5::OffsetHelper::order, ORDER_DCMP_IDX, Nektar::LibUtilities::ShapeTypeMap, Nektar::LibUtilities::SIZE_BasisType, Nektar::LibUtilities::SIZE_PointsType, and Nektar::LibUtilities::SIZE_ShapeType.

Referenced by v_Import().

◆ ImportHDF5FieldMetaData()

void Nektar::LibUtilities::FieldIOHdf5::ImportHDF5FieldMetaData ( DataSourceSharedPtr  dataSource,
FieldMetaDataMap fieldmetadatamap 
)
private

Import field metadata from dataSource.

Parameters
dataSourceInput datasource, which should be a H5DataSource.
fieldmetadatamapResulting field metadata from dataSource.

Definition at line 1410 of file FieldIOHdf5.cpp.

1412{
1414 std::static_pointer_cast<H5DataSource>(dataSource);
1415 H5::GroupSharedPtr master = hdf->Get()->OpenGroup("NEKTAR");
1416 H5::GroupSharedPtr metadata = master->OpenGroup("Metadata");
1417
1418 if (metadata)
1419 {
1420 H5::Group::AttrIterator param = metadata->attr_begin(),
1421 pEnd = metadata->attr_end();
1422 for (; param != pEnd; ++param)
1423 {
1424 std::string paramString = *param;
1425 if (paramString != "Provenance")
1426 {
1427 std::string paramBodyStr;
1428 metadata->GetAttribute(paramString, paramBodyStr);
1429 fieldmetadatamap[paramString] = paramBodyStr;
1430 }
1431 }
1432 }
1433}
std::shared_ptr< H5DataSource > H5DataSourceSharedPtr
Definition: FieldIOHdf5.h:86

Referenced by v_Import(), and v_ImportFieldMetaData().

◆ v_GetClassName()

const std::string & Nektar::LibUtilities::FieldIOHdf5::v_GetClassName ( ) const
overrideprotectedvirtual

Get class name.

Implements Nektar::LibUtilities::FieldIO.

Definition at line 1399 of file FieldIOHdf5.cpp.

1400{
1401 return className;
1402}
static std::string className
Name of class.
Definition: FieldIOHdf5.h:212

References className.

◆ v_Import()

void Nektar::LibUtilities::FieldIOHdf5::v_Import ( const std::string &  infilename,
std::vector< FieldDefinitionsSharedPtr > &  fielddefs,
std::vector< std::vector< NekDouble > > &  fielddata = NullVectorNekDoubleVector,
FieldMetaDataMap fieldinfomap = NullFieldMetaDataMap,
const Array< OneD, int > &  ElementIDs = NullInt1DArray 
)
overrideprotectedvirtual

Import a HDF5 format file.

Parameters
finfilenameInput filename
fielddefsField definitions of resulting field
fielddataField data of resulting field
fieldinfomapField metadata of resulting field
ElementIDsIf specified, contains the list of element IDs on this rank. The resulting field definitions will only contain data for the element IDs specified in this array.

Implements Nektar::LibUtilities::FieldIO.

Definition at line 885 of file FieldIOHdf5.cpp.

890{
891 std::stringstream prfx;
892 int nprocs = m_comm->GetSpaceComm()->GetSize();
893
894 // Set properties for parallel file access (if we're in parallel)
895 H5::PListSharedPtr parallelProps = H5::PList::Default();
898
899 if (nprocs > 1)
900 {
901 // Use MPI/O to access the file
902 parallelProps = H5::PList::FileAccess();
903 parallelProps->SetMpio(m_comm->GetSpaceComm());
904 // Use collective IO
905 readPL = H5::PList::DatasetXfer();
906 readPL->SetDxMpioCollective();
907 readPLInd = H5::PList::DatasetXfer();
908 readPLInd->SetDxMpioIndependent();
909 }
910
911 DataSourceSharedPtr dataSource =
912 H5DataSource::create(infilename, parallelProps);
913
914 // Open the root group of the hdf5 file
916 std::static_pointer_cast<H5DataSource>(dataSource);
917 ASSERTL1(h5, prfx.str() + "cannot open HDF5 file.");
918 H5::GroupSharedPtr root = h5->Get()->OpenGroup("NEKTAR");
919 ASSERTL1(root, prfx.str() + "cannot open root group.");
920
921 // Check format version
922 unsigned int formatVersion;
923 H5::Group::AttrIterator attrIt = root->attr_begin();
924 H5::Group::AttrIterator attrEnd = root->attr_end();
925 for (; attrIt != attrEnd; ++attrIt)
926 {
927 if (*attrIt == "FORMAT_VERSION")
928 {
929 break;
930 }
931 }
932
933 ASSERTL0(attrIt != attrEnd,
934 "Unable to determine Nektar++ HDF5 file version");
935 root->GetAttribute("FORMAT_VERSION", formatVersion);
936
937 ASSERTL0(formatVersion <= FORMAT_VERSION,
938 "File format if " + infilename +
939 " is higher than supported in "
940 "this version of Nektar++");
941
942 // Open the datasets
943 H5::DataSetSharedPtr decomps_dset = root->OpenDataSet("DECOMPOSITION");
944 ASSERTL1(decomps_dset, prfx.str() + "cannot open DECOMPOSITION dataset.");
945
946 H5::DataSetSharedPtr ids_dset = root->OpenDataSet("ELEMENTIDS");
947 ASSERTL1(ids_dset, prfx.str() + "cannot open ELEMENTIDS dataset.");
948
949 H5::DataSetSharedPtr data_dset = root->OpenDataSet("DATA");
950 ASSERTL1(data_dset, prfx.str() + "cannot open DATA dataset.");
951
952 // Open the dataset file spaces
953 H5::DataSpaceSharedPtr decomps_fspace = decomps_dset->GetSpace();
954 ASSERTL1(decomps_fspace,
955 prfx.str() + "cannot open DECOMPOSITION filespace.");
956
957 H5::DataSpaceSharedPtr ids_fspace = ids_dset->GetSpace();
958 ASSERTL1(ids_fspace, prfx.str() + "cannot open ELEMENTIDS filespace.");
959
960 H5::DataSpaceSharedPtr data_fspace = data_dset->GetSpace();
961 ASSERTL1(data_fspace, prfx.str() + "cannot open DATA filespace.");
962
963 // Read entire IDS data set to get list of global IDs.
964 std::vector<uint64_t> ids;
965
966 ids_dset->Read(ids, ids_fspace, readPL);
967
968 std::unordered_set<uint64_t> toread;
969 if (ElementIDs != NullInt1DArray)
970 {
971 for (uint64_t i = 0; i < ElementIDs.size(); ++i)
972 {
973 toread.insert(ElementIDs[i]);
974 }
975 }
976
977 std::vector<uint64_t> decomps;
978 decomps_dset->Read(decomps, decomps_fspace, readPL);
979
980 size_t nDecomps = decomps.size() / MAX_DCMPS;
981
982 // Mapping from each decomposition to offsets in the data array.
983 std::vector<OffsetHelper> decompsToOffsets(nDecomps);
984
985 // Mapping from each decomposition to a vector of element IDs. Note this has
986 // to be unsigned int, since that's what we use in FieldDefinitions.
987 std::map<uint64_t, std::vector<unsigned int>> decompsToElmts;
988
989 // Mapping from each group's hash to each of the decompositions.
990 std::map<uint64_t, std::set<uint64_t>> groupsToDecomps;
991
992 // True if we are pulling element IDs from ElementIDs.
993 bool selective = toread.size() > 0;
994
995 // Counters for data offsets
996 OffsetHelper running;
997
998 for (size_t i = 0, cnt = 0, cnt2 = 0; i < nDecomps; ++i, cnt += MAX_DCMPS)
999 {
1000 uint64_t nElmt = decomps[cnt + ELEM_DCMP_IDX];
1001 uint64_t groupHash = decomps[cnt + HASH_DCMP_IDX];
1002
1003 // Number of elements in this decomposition that this process needs
1004 uint64_t nElmtSelective = 0;
1005
1006 // Check if we should keep any elements in this decomposition
1007 if (selective)
1008 {
1009 for (size_t j = 0; j < nElmt; ++j)
1010 {
1011 uint64_t elmtId = ids[cnt2 + j];
1012 if (toread.find(elmtId) != toread.end())
1013 {
1014 nElmtSelective += 1;
1015 }
1016 }
1017 }
1018 else
1019 {
1020 nElmtSelective = nElmt;
1021 }
1022
1023 if (nElmtSelective > 0)
1024 {
1025 groupsToDecomps[groupHash].insert(i);
1026 }
1027
1028 // All element IDs in this decomposition
1029 std::vector<unsigned int> tmp2(nElmt);
1030 for (size_t j = 0; j < nElmt; ++j)
1031 {
1032 tmp2[j] = ids[cnt2 + j];
1033 }
1034
1035 cnt2 += nElmt;
1036
1037 decompsToElmts[i] = tmp2;
1038 decompsToOffsets[i] = running;
1039
1040 //
1041 // TODO: Will this update the OffsetHelper in decompsToOffsets? (copy
1042 // constructor?)
1043 //
1044 running.data += decomps[cnt + VAL_DCMP_IDX];
1045 running.order += decomps[cnt + ORDER_DCMP_IDX];
1046 running.homy += decomps[cnt + HOMY_DCMP_IDX];
1047 running.homz += decomps[cnt + HOMZ_DCMP_IDX];
1048 running.homs += decomps[cnt + HOMS_DCMP_IDX];
1049 }
1050
1051 for (auto &gIt : groupsToDecomps)
1052 {
1053 // Select region from dataset for this decomposition.
1054 for (auto &sIt : gIt.second)
1055 {
1056 // Convert group name to string
1057 auto groupName = std::to_string(gIt.first);
1058
1059 FieldDefinitionsSharedPtr fielddef =
1061 ImportFieldDef(readPLInd, root, decomps, sIt, decompsToOffsets[sIt],
1062 groupName, fielddef);
1063
1064 fielddef->m_elementIDs = decompsToElmts[sIt];
1065 //
1066 // TODO: Fix the case with 0 elements!
1067 //
1068 fielddefs.push_back(fielddef);
1069
1070 if (fielddata != NullVectorNekDoubleVector)
1071 {
1072 if (selective)
1073 {
1074 // Determine number of modes (coefficients) per element
1075 std::vector<unsigned int> coeffsPerElmt{
1076 GetNumberOfCoeffsPerElement(fielddef)};
1077
1078 // Selected element IDs
1079 std::vector<unsigned int> newElementIDs;
1080 std::vector<unsigned int> newNumModes;
1081 std::vector<hsize_t> dataIdxToRead;
1082
1083 {
1084 size_t offset = 0;
1085 size_t numbasis = fielddef->m_basis.size();
1086
1087 // Loop through all elements in this decomposition
1088
1089 for (size_t i = 0; i < fielddef->m_elementIDs.size();
1090 ++i)
1091 {
1092 // Check if we need data for this element
1093 if (toread.find(fielddef->m_elementIDs[i]) !=
1094 toread.end())
1095 {
1096 newElementIDs.push_back(
1097 fielddef->m_elementIDs[i]);
1098
1099 for (size_t j = 0; j < coeffsPerElmt[i]; ++j)
1100 {
1101 dataIdxToRead.push_back(
1102 decompsToOffsets[sIt].data + offset +
1103 j);
1104 }
1105
1106 if (fielddef->m_uniOrder == false)
1107 {
1108 for (size_t j = 0; j < numbasis; ++j)
1109 {
1110 newNumModes.push_back(
1111 fielddef
1112 ->m_numModes[i * numbasis + j]);
1113 }
1114 }
1115 }
1116
1117 offset += coeffsPerElmt[i];
1118 }
1119
1120 // Add indices for all remaining fields (variables)
1121 // We assume that all fields are stored with the same
1122 // polynomial order
1123 const size_t nDataPoints = dataIdxToRead.size();
1124
1125 for (size_t i = 1; i < fielddef->m_fields.size(); ++i)
1126 {
1127 for (size_t j = 0; j < nDataPoints; ++j)
1128 {
1129 dataIdxToRead.push_back(dataIdxToRead[j] +
1130 i * offset);
1131 }
1132 }
1133 }
1134
1135 fielddef->m_elementIDs = newElementIDs;
1136 if (fielddef->m_uniOrder == false)
1137 {
1138 fielddef->m_numModes = newNumModes;
1139 }
1140
1141 std::vector<NekDouble> decompFieldData;
1142
1143 data_fspace->ClearRange();
1144 data_fspace->SetSelection(dataIdxToRead.size(),
1145 dataIdxToRead);
1146
1147 data_dset->Read(decompFieldData, data_fspace, readPLInd);
1148 ASSERTL0(decompFieldData.size() ==
1149 CheckFieldDefinition(fielddef) *
1150 fielddef->m_fields.size(),
1151 "FieldIOHdf5: input data is not the same length "
1152 "as header information.");
1153 fielddata.push_back(decompFieldData);
1154 }
1155 else
1156 {
1157 std::vector<NekDouble> decompFieldData;
1158
1159 uint64_t nElemVals =
1160 decomps[sIt * MAX_DCMPS + VAL_DCMP_IDX];
1161 uint64_t nFieldVals = nElemVals;
1162
1163 data_fspace->ClearRange();
1164 data_fspace->SelectRange(decompsToOffsets[sIt].data,
1165 nFieldVals);
1166
1167 data_dset->Read(decompFieldData, data_fspace, readPLInd);
1168 ASSERTL0(decompFieldData.size() ==
1169 CheckFieldDefinition(fielddef) *
1170 fielddef->m_fields.size(),
1171 "FieldIOHdf5: input data is not the same length "
1172 "as header information.");
1173 fielddata.push_back(decompFieldData);
1174 }
1175 }
1176 }
1177 }
1178
1179 ImportHDF5FieldMetaData(dataSource, fieldinfomap);
1180 m_comm->GetSpaceComm()->Block();
1181}
static const unsigned int FORMAT_VERSION
Version of the Nektar++ HDF5 format, which is embedded into the main NEKTAR group as an attribute.
Definition: FieldIOHdf5.h:176
static const unsigned int HASH_DCMP_IDX
The hash of the field definition information, which defines the name of the attribute containing the ...
Definition: FieldIOHdf5.h:184
static const unsigned int ELEM_DCMP_IDX
A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of ...
Definition: FieldIOHdf5.h:178
void ImportFieldDef(H5::PListSharedPtr readPL, H5::GroupSharedPtr root, std::vector< uint64_t > &decomps, uint64_t decomp, OffsetHelper offset, std::string group, FieldDefinitionsSharedPtr def)
Import field definitions from a HDF5 document.
static const unsigned int VAL_DCMP_IDX
A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of ...
Definition: FieldIOHdf5.h:179
void ImportHDF5FieldMetaData(DataSourceSharedPtr dataSource, FieldMetaDataMap &fieldmetadatamap)
Import field metadata from dataSource.
int CheckFieldDefinition(const FieldDefinitionsSharedPtr &fielddefs)
Check field definitions for correctness and return storage size.
Definition: FieldIO.cpp:690
std::vector< unsigned int > GetNumberOfCoeffsPerElement(const FieldDefinitionsSharedPtr &fielddefs)
Compute number of data points needed to store expansion inside each element.
Definition: FieldIO.cpp:767
static PListSharedPtr DatasetXfer()
Properties for raw data transfer.
Definition: H5.cpp:109
static PListSharedPtr FileAccess()
Properties for file access.
Definition: H5.cpp:91
static PListSharedPtr Default()
Default options.
Definition: H5.cpp:73
static DataSourceSharedPtr create(const std::string &fn, H5::PListSharedPtr parallelProps)
Static constructor for this data source.
Definition: FieldIOHdf5.h:76
std::shared_ptr< PList > PListSharedPtr
Definition: H5.h:92
std::shared_ptr< DataSource > DataSourceSharedPtr
Definition: FieldIO.h:88
static std::vector< std::vector< NekDouble > > NullVectorNekDoubleVector
std::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
Definition: FieldIO.h:184
static Array< OneD, int > NullInt1DArray

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, ASSERTL1, Nektar::LibUtilities::FieldIO::CheckFieldDefinition(), Nektar::LibUtilities::H5DataSource::create(), Nektar::LibUtilities::FieldIOHdf5::OffsetHelper::data, Nektar::LibUtilities::H5::PList::DatasetXfer(), Nektar::LibUtilities::H5::PList::Default(), ELEM_DCMP_IDX, Nektar::LibUtilities::H5::PList::FileAccess(), FORMAT_VERSION, Nektar::LibUtilities::FieldIO::GetNumberOfCoeffsPerElement(), HASH_DCMP_IDX, Nektar::LibUtilities::FieldIOHdf5::OffsetHelper::homs, HOMS_DCMP_IDX, Nektar::LibUtilities::FieldIOHdf5::OffsetHelper::homy, HOMY_DCMP_IDX, Nektar::LibUtilities::FieldIOHdf5::OffsetHelper::homz, HOMZ_DCMP_IDX, ImportFieldDef(), ImportHDF5FieldMetaData(), Nektar::LibUtilities::FieldIO::m_comm, MAX_DCMPS, Nektar::NullInt1DArray, Nektar::LibUtilities::NullVectorNekDoubleVector, Nektar::LibUtilities::FieldIOHdf5::OffsetHelper::order, ORDER_DCMP_IDX, and VAL_DCMP_IDX.

◆ v_ImportFieldMetaData()

DataSourceSharedPtr Nektar::LibUtilities::FieldIOHdf5::v_ImportFieldMetaData ( const std::string &  filename,
FieldMetaDataMap fieldmetadatamap 
)
overrideprotectedvirtual

Import field metadata from filename and return the data source which wraps filename.

Parameters
filenameInput filename.
fieldmetadatamapResulting field metadata from dataSource.

Implements Nektar::LibUtilities::FieldIO.

Definition at line 1389 of file FieldIOHdf5.cpp.

1391{
1392 H5::PListSharedPtr parallelProps = H5::PList::Default();
1393 DataSourceSharedPtr ans = H5DataSource::create(filename, parallelProps);
1394 ImportHDF5FieldMetaData(ans, fieldmetadatamap);
1395 return ans;
1396}

References Nektar::LibUtilities::H5DataSource::create(), Nektar::LibUtilities::H5::PList::Default(), and ImportHDF5FieldMetaData().

◆ v_Write()

void Nektar::LibUtilities::FieldIOHdf5::v_Write ( const std::string &  outFile,
std::vector< FieldDefinitionsSharedPtr > &  fielddefs,
std::vector< std::vector< NekDouble > > &  fielddata,
const FieldMetaDataMap fieldmetadatamap = NullFieldMetaDataMap,
const bool  backup = false 
)
overrideprotectedvirtual

Write a HDF5 file to outFile given the field definitions fielddefs, field data fielddata and metadata fieldmetadatamap.

The writing strategy for HDF5 output is as follows:

  • Each rank determines the amount of data needed to be written into each dataset.
  • Each rank communicates its decomposition information to the root process.
  • The root processor initialises the output structure, writes the decomposition dataset and all the field definition information.
  • Other ranks may have field definitions that do not belong to the root process, in which case they open the file and append this (since attributes cannot be written in parallel).
  • Each of the other ranks writes their data contributions to the rest of the set.
Parameters
outFileOutput filename.
fielddefsInput field definitions.
fielddataInput field data.
fieldmetadatamapField metadata.

Implements Nektar::LibUtilities::FieldIO.

Definition at line 171 of file FieldIOHdf5.cpp.

176{
177 std::stringstream prfx;
178 prfx << m_comm->GetRank() << ": FieldIOHdf5::v_Write(): ";
179 double tm0 = 0.0, tm1 = 0.0;
180
181 if (m_comm->GetRank() == 0)
182 {
183 tm0 = m_comm->Wtime();
184 }
185
186 SetUpOutput(outFile, false, backup);
187
188 // We make a number of assumptions in this code:
189 // 1. All element ids have the same type: unsigned int
190 // 2. All elements within a given field have the same number of values
191 // 3. All element values have the same type, NekDouble
192
193 // Determine the root MPI process, i.e., the lowest ranked process handling
194 // nMaxFields fields, that will be responsible for writing our file.
195 ASSERTL1(fielddefs.size() == fielddata.size(),
196 prfx.str() + "fielddefs and fielddata have incompatible lengths.");
197
198 size_t nFields = fielddefs.size();
199 size_t nMaxFields = nFields;
200 m_comm->GetSpaceComm()->AllReduce(nMaxFields, LibUtilities::ReduceMax);
201
202 int root_rank = -1;
203 int nprocs = m_comm->GetSpaceComm()->GetSize();
204 bool amRoot = false;
205 LibUtilities::CommSharedPtr max_fields_comm;
206
207 if (nprocs > 1)
208 {
209 max_fields_comm = m_comm->GetSpaceComm()->CommCreateIf(
210 (nFields == nMaxFields) ? 1 : 0);
211 }
212 else
213 {
214 max_fields_comm = m_comm->GetSpaceComm();
215 }
216
217 if (max_fields_comm)
218 {
219 int rank = m_comm->GetSpaceComm()->GetRank();
220 root_rank = rank;
221 max_fields_comm->AllReduce(root_rank, LibUtilities::ReduceMin);
222 amRoot = (rank == root_rank);
223 if (!amRoot)
224 {
225 root_rank = -1;
226 }
227 }
228
229 m_comm->GetSpaceComm()->AllReduce(root_rank, LibUtilities::ReduceMax);
230 ASSERTL1(root_rank >= 0 && root_rank < nprocs,
231 prfx.str() + "invalid root rank.");
232
233 std::vector<uint64_t> decomps(nMaxFields * MAX_DCMPS, 0);
234 std::vector<uint64_t> all_hashes(nMaxFields * nprocs, 0);
235 std::vector<uint64_t> cnts(MAX_CNTS, 0);
236 std::vector<std::string> fieldNames(nFields);
237 std::vector<std::string> shapeStrings(nFields);
238 std::vector<std::vector<NekDouble>> homoLengths(nFields);
239 std::vector<std::vector<unsigned int>> homoSIDs(nFields), homoYIDs(nFields),
240 homoZIDs(nFields);
241 std::vector<std::vector<unsigned int>> numModesPerDirVar(nFields);
242 std::vector<std::string> numModesPerDirUni(nFields);
243
244 int homDim = -1;
245 int varOrder = 0;
246
247 for (int f = 0; f < nFields; ++f)
248 {
249 if (!fielddefs[f]->m_uniOrder)
250 {
251 varOrder = 1;
252 break;
253 }
254 }
255
256 m_comm->GetSpaceComm()->AllReduce(varOrder, LibUtilities::ReduceMax);
257
258 // Calculate the total number of elements handled by this MPI process and
259 // the total number of bytes required to store the elements. Base the name
260 // of each field on the hash of the field definition.
261 for (int f = 0; f < nFields; ++f)
262 {
263 ASSERTL1(fielddata[f].size() > 0,
264 prfx.str() +
265 "fielddata vector must contain at least one value.");
266 ASSERTL1(fielddata[f].size() == fielddefs[f]->m_fields.size() *
267 CheckFieldDefinition(fielddefs[f]),
268 prfx.str() + "fielddata vector has invalid size.");
269
270 std::size_t nFieldElems = fielddefs[f]->m_elementIDs.size();
271 std::size_t nElemVals = fielddata[f].size();
272
273 decomps[f * MAX_DCMPS + ELEM_DCMP_IDX] = nFieldElems;
274 decomps[f * MAX_DCMPS + VAL_DCMP_IDX] = nElemVals;
275
276 cnts[ELEM_CNT_IDX] += nFieldElems;
277 cnts[VAL_CNT_IDX] += nElemVals;
278
279 // Hash the field specification
280 std::stringstream hashStream;
281 size_t nSubFields = fielddefs[f]->m_fields.size();
282 for (size_t sf = 0; sf < nSubFields; ++sf)
283 {
284 hashStream << fielddefs[f]->m_fields[sf];
285 }
286
287 nSubFields = fielddefs[f]->m_basis.size();
288 for (size_t sf = 0; sf < nSubFields; ++sf)
289 {
290 hashStream << fielddefs[f]->m_basis[sf];
291 }
292
293 // Determine SHAPE attribute
294 std::stringstream shapeStringStream;
295 shapeStringStream << ShapeTypeMap[fielddefs[f]->m_shapeType];
296
297 if (fielddefs[f]->m_numHomogeneousDir > 0)
298 {
299 if (homDim == -1)
300 {
301 homDim = fielddefs[f]->m_numHomogeneousDir;
302 }
303
304 ASSERTL1(homDim == fielddefs[f]->m_numHomogeneousDir,
305 "HDF5 does not support variable homogeneous directions in "
306 "the same file.");
307
308 shapeStringStream << "-HomogenousExp"
309 << fielddefs[f]->m_numHomogeneousDir << "D";
310 }
311
312 if (fielddefs[f]->m_homoStrips)
313 {
314 shapeStringStream << "-Strips";
315 }
316
317 shapeStrings[f] = shapeStringStream.str();
318 hashStream << shapeStringStream.str();
319
320 // Determine HOMOGENEOUS attributes
321 if (fielddefs[f]->m_numHomogeneousDir)
322 {
323 nSubFields = fielddefs[f]->m_homogeneousLengths.size();
324 homoLengths[f].resize(nSubFields);
325 for (size_t sf = 0; sf < nSubFields; ++sf)
326 {
327 NekDouble len = fielddefs[f]->m_homogeneousLengths[sf];
328 hashStream << len;
329 homoLengths[f][sf] = len;
330 }
331
332 nSubFields = fielddefs[f]->m_homogeneousYIDs.size();
333 if (nSubFields > 0)
334 {
335 homoYIDs[f].resize(nSubFields);
336 decomps[f * MAX_DCMPS + HOMY_DCMP_IDX] = nSubFields;
337 cnts[HOMY_CNT_IDX] += nSubFields;
338 for (size_t sf = 0; sf < nSubFields; ++sf)
339 {
340 homoYIDs[f][sf] = fielddefs[f]->m_homogeneousYIDs[sf];
341 }
342 }
343
344 nSubFields = fielddefs[f]->m_homogeneousZIDs.size();
345 if (nSubFields > 0)
346 {
347 homoZIDs[f].resize(nSubFields);
348 decomps[f * MAX_DCMPS + HOMZ_DCMP_IDX] = nSubFields;
349 cnts[HOMZ_CNT_IDX] += nSubFields;
350 for (size_t sf = 0; sf < nSubFields; ++sf)
351 {
352 homoZIDs[f][sf] = fielddefs[f]->m_homogeneousZIDs[sf];
353 }
354 }
355
356 nSubFields = fielddefs[f]->m_homogeneousSIDs.size();
357 if (nSubFields > 0)
358 {
359 homoSIDs[f].resize(nSubFields);
360 decomps[f * MAX_DCMPS + HOMS_DCMP_IDX] = nSubFields;
361 cnts[HOMS_CNT_IDX] += nSubFields;
362 for (size_t sf = 0; sf < nSubFields; ++sf)
363 {
364 homoSIDs[f][sf] = fielddefs[f]->m_homogeneousSIDs[sf];
365 }
366 }
367 }
368
369 if (fielddefs[f]->m_uniOrder)
370 {
371 std::vector<unsigned int> elemModes(fielddefs[f]->m_basis.size());
372
373 for (size_t i = 0; i < fielddefs[f]->m_basis.size(); ++i)
374 {
375 elemModes[i] = fielddefs[f]->m_numModes[i];
376 }
377
378 if (varOrder)
379 {
380 for (size_t i = 0; i < nFieldElems; ++i)
381 {
382 std::copy(elemModes.begin(), elemModes.end(),
383 std::back_inserter(numModesPerDirVar[f]));
384 }
385 decomps[f * MAX_DCMPS + ORDER_DCMP_IDX] =
386 nFieldElems * elemModes.size();
387 cnts[ORDER_CNT_IDX] += nFieldElems * elemModes.size();
388 }
389 else
390 {
391 std::stringstream numModesStringStream;
392 numModesStringStream << "UNIORDER:";
393 for (size_t i = 0; i < elemModes.size(); i++)
394 {
395 if (i > 0)
396 {
397 numModesStringStream << ",";
398 }
399 numModesStringStream << elemModes[i];
400 }
401
402 numModesPerDirUni[f] = numModesStringStream.str();
403 hashStream << numModesPerDirUni[f];
404 }
405 }
406 else
407 {
408 numModesPerDirVar[f] = fielddefs[f]->m_numModes;
409 decomps[f * MAX_DCMPS + ORDER_DCMP_IDX] =
410 fielddefs[f]->m_numModes.size();
411 cnts[ORDER_CNT_IDX] += fielddefs[f]->m_numModes.size();
412 }
413
414 std::hash<std::string> string_hasher;
415 std::stringstream fieldNameStream;
416 uint64_t fieldDefHash = string_hasher(hashStream.str());
417
418 decomps[f * MAX_DCMPS + HASH_DCMP_IDX] = fieldDefHash;
419 all_hashes[m_comm->GetSpaceComm()->GetRank() * nMaxFields + f] =
420 fieldDefHash;
421
422 fieldNameStream << fieldDefHash;
423 fieldNames[f] = fieldNameStream.str();
424 }
425
426 // Gather information from all MPI processes
427 std::vector<uint64_t> all_cnts =
428 m_comm->GetSpaceComm()->Gather(root_rank, cnts);
429 std::vector<uint64_t> all_idxs(nprocs * MAX_IDXS, 0);
430 std::vector<uint64_t> all_decomps =
431 m_comm->GetSpaceComm()->Gather(root_rank, decomps);
432 std::vector<uint64_t> all_dsetsize(MAX_CNTS, 0);
433
434 // The root rank creates the file layout from scratch
435 if (amRoot)
436 {
437 H5::FileSharedPtr outfile = H5::File::Create(outFile, H5F_ACC_TRUNC);
438 ASSERTL1(outfile, prfx.str() + "cannot create HDF5 file.");
439 H5::GroupSharedPtr root = outfile->CreateGroup("NEKTAR");
440 ASSERTL1(root, prfx.str() + "cannot create root group.");
441 TagWriterSharedPtr info_writer(new H5TagWriter(root));
442 AddInfoTag(info_writer, fieldmetadatamap);
443
444 // Record file format version as attribute in main group.
445 root->SetAttribute("FORMAT_VERSION", FORMAT_VERSION);
446
447 // Calculate the indexes to be used by each MPI process when reading the
448 // IDS and DATA datasets
449 size_t nTotElems = 0, nTotVals = 0, nTotOrder = 0;
450 size_t nTotHomY = 0, nTotHomZ = 0, nTotHomS = 0;
451 for (int r = 0; r < nprocs; ++r)
452 {
453 all_idxs[r * MAX_IDXS + IDS_IDX_IDX] = nTotElems;
454 all_idxs[r * MAX_IDXS + DATA_IDX_IDX] = nTotVals;
455 all_idxs[r * MAX_IDXS + ORDER_IDX_IDX] = nTotOrder;
456 all_idxs[r * MAX_IDXS + HOMY_IDX_IDX] = nTotHomY;
457 all_idxs[r * MAX_IDXS + HOMZ_IDX_IDX] = nTotHomZ;
458 all_idxs[r * MAX_IDXS + HOMS_IDX_IDX] = nTotHomS;
459
460 nTotElems += all_cnts[r * MAX_CNTS + ELEM_CNT_IDX];
461 nTotVals += all_cnts[r * MAX_CNTS + VAL_CNT_IDX];
462 nTotOrder += all_cnts[r * MAX_CNTS + ORDER_CNT_IDX];
463 nTotHomY += all_cnts[r * MAX_CNTS + HOMY_CNT_IDX];
464 nTotHomZ += all_cnts[r * MAX_CNTS + HOMZ_CNT_IDX];
465 nTotHomS += all_cnts[r * MAX_CNTS + HOMS_CNT_IDX];
466 }
467
468 all_dsetsize[ELEM_CNT_IDX] = nTotElems;
469 all_dsetsize[VAL_CNT_IDX] = nTotVals;
470 all_dsetsize[ORDER_CNT_IDX] = nTotOrder;
471 all_dsetsize[HOMY_CNT_IDX] = nTotHomY;
472 all_dsetsize[HOMZ_CNT_IDX] = nTotHomZ;
473 all_dsetsize[HOMS_CNT_IDX] = nTotHomS;
474
475 // Create DECOMPOSITION dataset: basic field info for each MPI process
476 H5::DataTypeSharedPtr decomps_type =
477 H5::DataType::OfObject(all_decomps[0]);
478 H5::DataSpaceSharedPtr decomps_space =
479 H5::DataSpace::OneD(all_decomps.size());
480 H5::DataSetSharedPtr decomps_dset =
481 root->CreateDataSet("DECOMPOSITION", decomps_type, decomps_space);
482 ASSERTL1(decomps_dset,
483 prfx.str() + "cannot create DECOMPOSITION dataset.");
484
485 // Create IDS dataset: element ids
486 H5::DataTypeSharedPtr ids_type =
487 H5::DataType::OfObject(fielddefs[0]->m_elementIDs[0]);
488 H5::DataSpaceSharedPtr ids_space = H5::DataSpace::OneD(nTotElems);
489 H5::DataSetSharedPtr ids_dset =
490 root->CreateDataSet("ELEMENTIDS", ids_type, ids_space);
491 ASSERTL1(ids_dset, prfx.str() + "cannot create ELEMENTIDS dataset.");
492
493 // Create DATA dataset: element data
494 H5::DataTypeSharedPtr data_type =
495 H5::DataType::OfObject(fielddata[0][0]);
496 H5::DataSpaceSharedPtr data_space = H5::DataSpace::OneD(nTotVals);
497 H5::DataSetSharedPtr data_dset =
498 root->CreateDataSet("DATA", data_type, data_space);
499 ASSERTL1(data_dset, prfx.str() + "cannot create DATA dataset.");
500
501 // Create HOMOGENEOUSYIDS dataset: homogeneous y-plane IDs
502 if (nTotHomY > 0)
503 {
504 H5::DataTypeSharedPtr homy_type =
505 H5::DataType::OfObject(homoYIDs[0][0]);
506 H5::DataSpaceSharedPtr homy_space = H5::DataSpace::OneD(nTotHomY);
507 H5::DataSetSharedPtr homy_dset =
508 root->CreateDataSet("HOMOGENEOUSYIDS", homy_type, homy_space);
509 ASSERTL1(homy_dset,
510 prfx.str() + "cannot create HOMOGENEOUSYIDS dataset.");
511 }
512
513 // Create HOMOGENEOUSYIDS dataset: homogeneous z-plane IDs
514 if (nTotHomZ > 0)
515 {
516 H5::DataTypeSharedPtr homz_type =
517 H5::DataType::OfObject(homoZIDs[0][0]);
518 H5::DataSpaceSharedPtr homz_space = H5::DataSpace::OneD(nTotHomZ);
519 H5::DataSetSharedPtr homz_dset =
520 root->CreateDataSet("HOMOGENEOUSZIDS", homz_type, homz_space);
521 ASSERTL1(homz_dset,
522 prfx.str() + "cannot create HOMOGENEOUSZIDS dataset.");
523 }
524
525 // Create HOMOGENEOUSSIDS dataset: homogeneous strip IDs
526 if (nTotHomS > 0)
527 {
528 H5::DataTypeSharedPtr homs_type =
529 H5::DataType::OfObject(homoSIDs[0][0]);
530 H5::DataSpaceSharedPtr homs_space = H5::DataSpace::OneD(nTotHomS);
531 H5::DataSetSharedPtr homs_dset =
532 root->CreateDataSet("HOMOGENEOUSSIDS", homs_type, homs_space);
533 ASSERTL1(homs_dset,
534 prfx.str() + "cannot create HOMOGENEOUSSIDS dataset.");
535 }
536
537 // Create POLYORDERS dataset: elemental polynomial orders
538 if (varOrder)
539 {
540 H5::DataTypeSharedPtr order_type =
541 H5::DataType::OfObject(numModesPerDirVar[0][0]);
542 H5::DataSpaceSharedPtr order_space = H5::DataSpace::OneD(nTotOrder);
543 H5::DataSetSharedPtr order_dset =
544 root->CreateDataSet("POLYORDERS", order_type, order_space);
545 ASSERTL1(order_dset,
546 prfx.str() + "cannot create POLYORDERS dataset.");
547 }
548 }
549
550 m_comm->GetSpaceComm()->Bcast(all_dsetsize, root_rank);
551
552 // Datasets, root group and HDF5 file are all closed automatically since
553 // they are now out of scope. Now we need to determine which process will
554 // write the group representing the field description in the HDF5 file. This
555 // next block of code performs this by finding all unique hashes and then
556 // determining one process that will create (possibly more than one) group
557 // for that hash. An alternative would be to communicate the field
558 // information to the root processor, but this is a bit convoluted.
559
560 // This set stores the unique hashes.
561 std::set<uint64_t> hashToProc;
562
563 // This map takes ranks to hashes this process will write.
564 std::map<int, std::vector<uint64_t>> writingProcs;
565
566 // Gather all field hashes to every processor.
567 m_comm->GetSpaceComm()->AllReduce(all_hashes, LibUtilities::ReduceMax);
568
569 for (int n = 0; n < nprocs; ++n)
570 {
571 for (int i = 0; i < nMaxFields; ++i)
572 {
573 uint64_t hash = all_hashes[n * nMaxFields + i];
574
575 // Note hash can be zero if, on this process, nFields < nMaxFields.
576 if (hashToProc.find(hash) != hashToProc.end() || hash == 0)
577 {
578 continue;
579 }
580 hashToProc.insert(hash);
581 writingProcs[n].push_back(hash);
582 }
583 }
584
585 // Having constructed the map, go ahead and write the attributes out.
586 for (auto &sIt : writingProcs)
587 {
588 int rank = sIt.first;
589
590 // Write out this rank's groups.
591 if (m_comm->GetSpaceComm()->GetRank() == rank)
592 {
595
596 // Reopen the file
597 H5::FileSharedPtr outfile =
598 H5::File::Open(outFile, H5F_ACC_RDWR, serialProps);
599 ASSERTL1(outfile, prfx.str() + "cannot open HDF5 file.");
600 H5::GroupSharedPtr root = outfile->OpenGroup("NEKTAR");
601 ASSERTL1(root, prfx.str() + "cannot open root group.");
602
603 // Write a HDF5 group for each field
604 hashToProc.clear();
605 for (int i = 0; i < sIt.second.size(); ++i)
606 {
607 for (int f = 0; f < nFields; ++f)
608 {
609 if (sIt.second[i] !=
610 all_hashes[m_comm->GetSpaceComm()->GetRank() *
611 nMaxFields +
612 f] ||
613 hashToProc.find(sIt.second[i]) != hashToProc.end())
614 {
615 continue;
616 }
617
618 hashToProc.insert(sIt.second[i]);
619
620 // Just in case we've already written this
621 H5::GroupSharedPtr field_group =
622 root->CreateGroup(fieldNames[f]);
623 ASSERTL1(field_group,
624 prfx.str() + "cannot create field group.");
625 field_group->SetAttribute("FIELDS", fielddefs[f]->m_fields);
626 field_group->SetAttribute("BASIS", fielddefs[f]->m_basis);
627 field_group->SetAttribute("SHAPE", shapeStrings[f]);
628
629 if (homoLengths[f].size() > 0)
630 {
631 field_group->SetAttribute("HOMOGENEOUSLENGTHS",
632 homoLengths[f]);
633 }
634
635 // If the field has only uniform order, we write the order
636 // into the NUMMODESPERDIR attribute. Otherwise, we'll go
637 // ahead and assume everything is mixed and fix this in the
638 // read later if required.
639 if (!varOrder)
640 {
641 field_group->SetAttribute("NUMMODESPERDIR",
642 numModesPerDirUni[f]);
643 }
644 else
645 {
646 std::string numModesPerDir = "MIXORDER";
647 field_group->SetAttribute("NUMMODESPERDIR",
648 numModesPerDir);
649 }
650 }
651 }
652 }
653
654 // We block to avoid more than one processor opening the file at a time.
655 m_comm->GetSpaceComm()->Block();
656 }
657
658 // Write the DECOMPOSITION dataset
659 if (amRoot)
660 {
663
664 // Reopen the file
665 H5::FileSharedPtr outfile =
666 H5::File::Open(outFile, H5F_ACC_RDWR, serialProps);
667 ASSERTL1(outfile, prfx.str() + "cannot open HDF5 file.");
668 H5::GroupSharedPtr root = outfile->OpenGroup("NEKTAR");
669 ASSERTL1(root, prfx.str() + "cannot open root group.");
670
671 // Write the DECOMPOSITION dataset
672 H5::DataSetSharedPtr decomps_dset = root->OpenDataSet("DECOMPOSITION");
673 ASSERTL1(decomps_dset,
674 prfx.str() + "cannot open DECOMPOSITION dataset.");
675
676 H5::DataSpaceSharedPtr decomps_fspace = decomps_dset->GetSpace();
677 ASSERTL1(decomps_fspace,
678 prfx.str() + "cannot open DECOMPOSITION filespace.");
679
680 decomps_fspace->SelectRange(0, all_decomps.size());
681 decomps_dset->Write(all_decomps, decomps_fspace, writeSR);
682 }
683
684 // Initialise the dataset indexes for all MPI processes
685 std::vector<uint64_t> idx =
686 m_comm->GetSpaceComm()->Scatter(root_rank, all_idxs);
687 uint64_t ids_i = idx[IDS_IDX_IDX];
688 uint64_t data_i = idx[DATA_IDX_IDX];
689 uint64_t order_i = idx[ORDER_IDX_IDX];
690 uint64_t homy_i = idx[HOMY_IDX_IDX];
691 uint64_t homz_i = idx[HOMZ_IDX_IDX];
692 uint64_t homs_i = idx[HOMS_IDX_IDX];
693
694 // Set properties for parallel file access (if we're in parallel)
695 H5::PListSharedPtr parallelProps = H5::PList::Default();
697 if (nprocs > 1)
698 {
699 // Use MPI/O to access the file
700 parallelProps = H5::PList::FileAccess();
701 parallelProps->SetMpio(m_comm->GetSpaceComm());
702 // Use collective IO
703 writePL = H5::PList::DatasetXfer();
704 writePL->SetDxMpioCollective();
705 }
706
707 // Reopen the file
708 H5::FileSharedPtr outfile =
709 H5::File::Open(outFile, H5F_ACC_RDWR, parallelProps);
710 ASSERTL1(outfile, prfx.str() + "cannot open HDF5 file.");
711 H5::GroupSharedPtr root = outfile->OpenGroup("NEKTAR");
712 ASSERTL1(root, prfx.str() + "cannot open root group.");
713
714 m_comm->GetSpaceComm()->Block();
715
716 // all HDF5 groups have now been created. Open the IDS dataset and
717 // associated data space
718 H5::DataSetSharedPtr ids_dset = root->OpenDataSet("ELEMENTIDS");
719 ASSERTL1(ids_dset, prfx.str() + "cannot open ELEMENTIDS dataset.");
720 H5::DataSpaceSharedPtr ids_fspace = ids_dset->GetSpace();
721 ASSERTL1(ids_fspace, prfx.str() + "cannot open ELEMENTIDS filespace.");
722
723 // Open the DATA dataset and associated data space
724 H5::DataSetSharedPtr data_dset = root->OpenDataSet("DATA");
725 ASSERTL1(data_dset, prfx.str() + "cannot open DATA dataset.");
726 H5::DataSpaceSharedPtr data_fspace = data_dset->GetSpace();
727 ASSERTL1(data_fspace, prfx.str() + "cannot open DATA filespace.");
728
729 // Open the optional datasets and data spaces.
730 H5::DataSetSharedPtr order_dset, homy_dset, homz_dset, homs_dset;
731 H5::DataSpaceSharedPtr order_fspace, homy_fspace, homz_fspace, homs_fspace;
732
733 if (all_dsetsize[ORDER_CNT_IDX])
734 {
735 order_dset = root->OpenDataSet("POLYORDERS");
736 ASSERTL1(order_dset, prfx.str() + "cannot open POLYORDERS dataset.");
737 order_fspace = order_dset->GetSpace();
738 ASSERTL1(order_fspace,
739 prfx.str() + "cannot open POLYORDERS filespace.");
740 }
741
742 if (all_dsetsize[HOMY_CNT_IDX])
743 {
744 homy_dset = root->OpenDataSet("HOMOGENEOUSYIDS");
745 ASSERTL1(homy_dset,
746 prfx.str() + "cannot open HOMOGENEOUSYIDS dataset.");
747 homy_fspace = homy_dset->GetSpace();
748 ASSERTL1(homy_fspace,
749 prfx.str() + "cannot open HOMOGENEOUSYIDS filespace.");
750 }
751
752 if (all_dsetsize[HOMZ_CNT_IDX])
753 {
754 homz_dset = root->OpenDataSet("HOMOGENEOUSZIDS");
755 ASSERTL1(homz_dset,
756 prfx.str() + "cannot open HOMOGENEOUSZIDS dataset.");
757 homz_fspace = homz_dset->GetSpace();
758 ASSERTL1(homz_fspace,
759 prfx.str() + "cannot open HOMOGENEOUSZIDS filespace.");
760 }
761
762 if (all_dsetsize[HOMS_CNT_IDX])
763 {
764 homs_dset = root->OpenDataSet("HOMOGENEOUSSIDS");
765 ASSERTL1(homs_dset,
766 prfx.str() + "cannot open HOMOGENEOUSSIDS dataset.");
767 homs_fspace = homs_dset->GetSpace();
768 ASSERTL1(homs_fspace,
769 prfx.str() + "cannot open HOMOGENEOUSSIDS filespace.");
770 }
771
772 // Write the data
773 for (int f = 0; f < nFields; ++f)
774 {
775 // write the element ids
776 std::size_t nFieldElems = fielddefs[f]->m_elementIDs.size();
777 ids_fspace->SelectRange(ids_i, nFieldElems);
778 ids_dset->Write(fielddefs[f]->m_elementIDs, ids_fspace, writePL);
779 ids_i += nFieldElems;
780
781 // write the element values
782 std::size_t nFieldVals = fielddata[f].size();
783 data_fspace->SelectRange(data_i, nFieldVals);
784 data_dset->Write(fielddata[f], data_fspace, writePL);
785 data_i += nFieldVals;
786 }
787
788 if (order_dset)
789 {
790 for (int f = 0; f < nFields; ++f)
791 {
792 std::size_t nOrders = numModesPerDirVar[f].size();
793 order_fspace->SelectRange(order_i, nOrders);
794 order_dset->Write(numModesPerDirVar[f], order_fspace, writePL);
795 order_i += nOrders;
796 }
797 }
798
799 if (homy_dset)
800 {
801 for (int f = 0; f < nFields; ++f)
802 {
803 std::size_t nYIDs = homoYIDs[f].size();
804 homy_fspace->SelectRange(homy_i, nYIDs);
805 homy_dset->Write(homoYIDs[f], homy_fspace, writePL);
806 homy_i += nYIDs;
807 }
808 }
809
810 if (homz_dset)
811 {
812 for (int f = 0; f < nFields; ++f)
813 {
814 std::size_t nZIDs = homoZIDs[f].size();
815 homz_fspace->SelectRange(homz_i, nZIDs);
816 homz_dset->Write(homoZIDs[f], homz_fspace, writePL);
817 homz_i += nZIDs;
818 }
819 }
820
821 if (homs_dset)
822 {
823 for (int f = 0; f < nFields; ++f)
824 {
825 std::size_t nSIDs = homoSIDs[f].size();
826 homs_fspace->SelectRange(homs_i, nSIDs);
827 homs_dset->Write(homoSIDs[f], homs_fspace, writePL);
828 homs_i += nSIDs;
829 }
830 }
831
832 for (int f = nFields; f < nMaxFields; ++f)
833 {
834 // this MPI process is handling fewer than nMaxFields fields
835 // so, since this is a collective operation
836 // just rewrite the element ids and values of the last field
837 ids_dset->Write(fielddefs[nFields - 1]->m_elementIDs, ids_fspace,
838 writePL);
839 data_dset->Write(fielddata[nFields - 1], data_fspace, writePL);
840
841 if (order_dset)
842 {
843 order_dset->Write(numModesPerDirVar[nFields - 1], order_fspace,
844 writePL);
845 }
846
847 if (homy_dset)
848 {
849 homy_dset->Write(homoYIDs[nFields - 1], homy_fspace, writePL);
850 }
851
852 if (homz_dset)
853 {
854 homz_dset->Write(homoZIDs[nFields - 1], homz_fspace, writePL);
855 }
856
857 if (homs_dset)
858 {
859 homs_dset->Write(homoSIDs[nFields - 1], homs_fspace, writePL);
860 }
861 }
862
863 m_comm->GetSpaceComm()->Block();
864
865 // all data has been written
866 if (m_comm->GetRank() == 0)
867 {
868 tm1 = m_comm->Wtime();
869 std::cout << " (" << tm1 - tm0 << "s, HDF5)" << std::endl;
870 }
871}
static const unsigned int ELEM_CNT_IDX
A helper for FieldIOHdf5::v_Write. Describes the position of the number of elements in the cnt array.
Definition: FieldIOHdf5.h:187
static const unsigned int MAX_CNTS
A helper for FieldIOHdf5::v_Write. Describes the maximum number of items in the cnt array per field d...
Definition: FieldIOHdf5.h:193
static const unsigned int HOMY_CNT_IDX
A helper for FieldIOHdf5::v_Write. Describes the position of the number of homogeneous y-planes in th...
Definition: FieldIOHdf5.h:190
static const unsigned int HOMZ_CNT_IDX
A helper for FieldIOHdf5::v_Write. Describes the position of the number of homogeneous z-planes in th...
Definition: FieldIOHdf5.h:191
static const unsigned int HOMS_CNT_IDX
A helper for FieldIOHdf5::v_Write. Describes the position of the number of homogeneous strips in the ...
Definition: FieldIOHdf5.h:192
static const unsigned int VAL_CNT_IDX
A helper for FieldIOHdf5::v_Write. Describes the position of the number of data points in the cnt arr...
Definition: FieldIOHdf5.h:188
static const unsigned int IDS_IDX_IDX
A helper for FieldIOHdf5::v_Write. Describes the position of the element IDs within the indexing set.
Definition: FieldIOHdf5.h:195
static const unsigned int ORDER_CNT_IDX
A helper for FieldIOHdf5::v_Write. Describes the position of the number of order points in the cnt ar...
Definition: FieldIOHdf5.h:189
static const unsigned int ORDER_IDX_IDX
A helper for FieldIOHdf5::v_Write. Describes the position of the element order within the indexing se...
Definition: FieldIOHdf5.h:197
static const unsigned int MAX_IDXS
A helper for FieldIOHdf5::v_Write. Describes the maximum number of items in the indexing set.
Definition: FieldIOHdf5.h:201
static const unsigned int DATA_IDX_IDX
A helper for FieldIOHdf5::v_Write. Describes the position of the data size within the indexing set.
Definition: FieldIOHdf5.h:196
static const unsigned int HOMS_IDX_IDX
A helper for FieldIOHdf5::v_Write. Describes the position of the number of homogeneous strips within ...
Definition: FieldIOHdf5.h:200
static const unsigned int HOMZ_IDX_IDX
A helper for FieldIOHdf5::v_Write. Describes the position of the number of z-planes within the indexi...
Definition: FieldIOHdf5.h:199
static const unsigned int HOMY_IDX_IDX
A helper for FieldIOHdf5::v_Write. Describes the position of the number of y-planes within the indexi...
Definition: FieldIOHdf5.h:198
std::string SetUpOutput(const std::string outname, bool perRank, bool backup=false)
Set up the filesystem ready for output.
Definition: FieldIO.cpp:404
static void AddInfoTag(TagWriterSharedPtr root, const FieldMetaDataMap &fieldmetadatamap)
Add provenance information to the field metadata map.
Definition: FieldIO.cpp:342
static DataSpaceSharedPtr OneD(hsize_t size)
Definition: H5.cpp:411
static DataTypeSharedPtr OfObject(const T &obj)
Definition: H5.h:397
static FileSharedPtr Open(const std::string &filename, unsigned mode, PListSharedPtr accessPL=PList::Default())
Definition: H5.cpp:636
static FileSharedPtr Create(const std::string &filename, unsigned mode, PListSharedPtr createPL=PList::Default(), PListSharedPtr accessPL=PList::Default())
Definition: H5.cpp:627
def copy(self)
Definition: pycml.py:2663
std::shared_ptr< File > FileSharedPtr
Definition: H5.h:88
std::shared_ptr< DataType > DataTypeSharedPtr
Definition: H5.h:74
std::shared_ptr< TagWriter > TagWriterSharedPtr
Definition: FieldIO.h:78
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
double NekDouble

References Nektar::LibUtilities::FieldIO::AddInfoTag(), ASSERTL1, Nektar::LibUtilities::FieldIO::CheckFieldDefinition(), CellMLToNektar.pycml::copy(), Nektar::LibUtilities::H5::File::Create(), DATA_IDX_IDX, Nektar::LibUtilities::H5::PList::DatasetXfer(), Nektar::LibUtilities::H5::PList::Default(), ELEM_CNT_IDX, ELEM_DCMP_IDX, Nektar::LibUtilities::H5::PList::FileAccess(), FORMAT_VERSION, HASH_DCMP_IDX, HOMS_CNT_IDX, HOMS_DCMP_IDX, HOMS_IDX_IDX, HOMY_CNT_IDX, HOMY_DCMP_IDX, HOMY_IDX_IDX, HOMZ_CNT_IDX, HOMZ_DCMP_IDX, HOMZ_IDX_IDX, IDS_IDX_IDX, Nektar::LibUtilities::FieldIO::m_comm, MAX_CNTS, MAX_DCMPS, MAX_IDXS, Nektar::LibUtilities::H5::DataType::OfObject(), Nektar::LibUtilities::H5::DataSpace::OneD(), Nektar::LibUtilities::H5::File::Open(), ORDER_CNT_IDX, ORDER_DCMP_IDX, ORDER_IDX_IDX, Nektar::LibUtilities::ReduceMax, Nektar::LibUtilities::ReduceMin, Nektar::LibUtilities::FieldIO::SetUpOutput(), Nektar::LibUtilities::ShapeTypeMap, VAL_CNT_IDX, and VAL_DCMP_IDX.

Member Data Documentation

◆ className

std::string Nektar::LibUtilities::FieldIOHdf5::className
static
Initial value:
=
"Hdf5", FieldIOHdf5::create, "HDF5-based output of field data.")
static FieldIOSharedPtr create(LibUtilities::CommSharedPtr pComm, bool sharedFilesystem)
Creates an instance of this class.
Definition: FieldIOHdf5.h:204
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
FieldIOFactory & GetFieldIOFactory()
Returns the FieldIO factory.
Definition: FieldIO.cpp:70

Name of class.

Definition at line 212 of file FieldIOHdf5.h.

Referenced by v_GetClassName().

◆ DATA_IDX_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::DATA_IDX_IDX = 1
static

A helper for FieldIOHdf5::v_Write. Describes the position of the data size within the indexing set.

Definition at line 196 of file FieldIOHdf5.h.

Referenced by v_Write().

◆ ELEM_CNT_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::ELEM_CNT_IDX = 0
static

A helper for FieldIOHdf5::v_Write. Describes the position of the number of elements in the cnt array.

Definition at line 187 of file FieldIOHdf5.h.

Referenced by v_Write().

◆ ELEM_DCMP_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::ELEM_DCMP_IDX = 0
static

A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of elements in decomposition (i.e. field definition).

Definition at line 178 of file FieldIOHdf5.h.

Referenced by v_Import(), and v_Write().

◆ FORMAT_VERSION

const unsigned int Nektar::LibUtilities::FieldIOHdf5::FORMAT_VERSION = 1
static

Version of the Nektar++ HDF5 format, which is embedded into the main NEKTAR group as an attribute.

Definition at line 176 of file FieldIOHdf5.h.

Referenced by v_Import(), and v_Write().

◆ HASH_DCMP_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::HASH_DCMP_IDX = 6
static

The hash of the field definition information, which defines the name of the attribute containing the field definition itself.

Definition at line 184 of file FieldIOHdf5.h.

Referenced by v_Import(), and v_Write().

◆ HOMS_CNT_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::HOMS_CNT_IDX = 5
static

A helper for FieldIOHdf5::v_Write. Describes the position of the number of homogeneous strips in the cnt array.

Definition at line 192 of file FieldIOHdf5.h.

Referenced by v_Write().

◆ HOMS_DCMP_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::HOMS_DCMP_IDX = 5
static

A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of the number of strips for homogeneous simulations.

Definition at line 183 of file FieldIOHdf5.h.

Referenced by ImportFieldDef(), v_Import(), and v_Write().

◆ HOMS_IDX_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::HOMS_IDX_IDX = 5
static

A helper for FieldIOHdf5::v_Write. Describes the position of the number of homogeneous strips within the indexing set.

Definition at line 200 of file FieldIOHdf5.h.

Referenced by v_Write().

◆ HOMY_CNT_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::HOMY_CNT_IDX = 3
static

A helper for FieldIOHdf5::v_Write. Describes the position of the number of homogeneous y-planes in the cnt array.

Definition at line 190 of file FieldIOHdf5.h.

Referenced by v_Write().

◆ HOMY_DCMP_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::HOMY_DCMP_IDX = 3
static

A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of the number of y-planes for homogeneous simulations.

Definition at line 181 of file FieldIOHdf5.h.

Referenced by ImportFieldDef(), v_Import(), and v_Write().

◆ HOMY_IDX_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::HOMY_IDX_IDX = 3
static

A helper for FieldIOHdf5::v_Write. Describes the position of the number of y-planes within the indexing set.

Definition at line 198 of file FieldIOHdf5.h.

Referenced by v_Write().

◆ HOMZ_CNT_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::HOMZ_CNT_IDX = 4
static

A helper for FieldIOHdf5::v_Write. Describes the position of the number of homogeneous z-planes in the cnt array.

Definition at line 191 of file FieldIOHdf5.h.

Referenced by v_Write().

◆ HOMZ_DCMP_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::HOMZ_DCMP_IDX = 4
static

A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of the number of z-planes for homogeneous simulations.

Definition at line 182 of file FieldIOHdf5.h.

Referenced by ImportFieldDef(), v_Import(), and v_Write().

◆ HOMZ_IDX_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::HOMZ_IDX_IDX = 4
static

A helper for FieldIOHdf5::v_Write. Describes the position of the number of z-planes within the indexing set.

Definition at line 199 of file FieldIOHdf5.h.

Referenced by v_Write().

◆ IDS_IDX_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::IDS_IDX_IDX = 0
static

A helper for FieldIOHdf5::v_Write. Describes the position of the element IDs within the indexing set.

Definition at line 195 of file FieldIOHdf5.h.

Referenced by v_Write().

◆ MAX_CNTS

const unsigned int Nektar::LibUtilities::FieldIOHdf5::MAX_CNTS = FieldIOHdf5::HOMS_CNT_IDX + 1
static

A helper for FieldIOHdf5::v_Write. Describes the maximum number of items in the cnt array per field definition.

Definition at line 193 of file FieldIOHdf5.h.

Referenced by v_Write().

◆ MAX_DCMPS

const unsigned int Nektar::LibUtilities::FieldIOHdf5::MAX_DCMPS = FieldIOHdf5::HASH_DCMP_IDX + 1
static

A helper for FieldIOHdf5::v_Write. Describes the maximum number of items in the decomposition per field definition.

Definition at line 185 of file FieldIOHdf5.h.

Referenced by ImportFieldDef(), v_Import(), and v_Write().

◆ MAX_IDXS

const unsigned int Nektar::LibUtilities::FieldIOHdf5::MAX_IDXS = FieldIOHdf5::HOMS_IDX_IDX + 1
static

A helper for FieldIOHdf5::v_Write. Describes the maximum number of items in the indexing set.

Definition at line 201 of file FieldIOHdf5.h.

Referenced by v_Write().

◆ ORDER_CNT_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::ORDER_CNT_IDX = 2
static

A helper for FieldIOHdf5::v_Write. Describes the position of the number of order points in the cnt array.

Definition at line 189 of file FieldIOHdf5.h.

Referenced by v_Write().

◆ ORDER_DCMP_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::ORDER_DCMP_IDX = 2
static

A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of elements multiplied by the dimension of the element, giving number of modes when variable polynomial order is defined.

Definition at line 180 of file FieldIOHdf5.h.

Referenced by ImportFieldDef(), v_Import(), and v_Write().

◆ ORDER_IDX_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::ORDER_IDX_IDX = 2
static

A helper for FieldIOHdf5::v_Write. Describes the position of the element order within the indexing set.

Definition at line 197 of file FieldIOHdf5.h.

Referenced by v_Write().

◆ VAL_CNT_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::VAL_CNT_IDX = 1
static

A helper for FieldIOHdf5::v_Write. Describes the position of the number of data points in the cnt array.

Definition at line 188 of file FieldIOHdf5.h.

Referenced by v_Write().

◆ VAL_DCMP_IDX

const unsigned int Nektar::LibUtilities::FieldIOHdf5::VAL_DCMP_IDX = 1
static

A helper for FieldIOHdf5::v_Write and FieldIOHdf5::v_Import. Describes the position of the number of data points in decomposition (i.e. field definition).

Definition at line 179 of file FieldIOHdf5.h.

Referenced by v_Import(), and v_Write().