Nektar++
MeshGraphHDF5.cpp
Go to the documentation of this file.
1 ////////////////////////////////////////////////////////////////////////////////
2 //
3 // File: MeshGraphHDF5.cpp
4 //
5 // For more information, please see: http://www.nektar.info/
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: HDF5-based mesh format for Nektar++.
32 //
33 ////////////////////////////////////////////////////////////////////////////////
34 
35 #include <type_traits>
36 #include <tinyxml.h>
37 #include <boost/core/ignore_unused.hpp>
38 #include <boost/algorithm/string.hpp>
39 #include <boost/filesystem.hpp>
40 
45 
46 #define TIME_RESULT(verb, msg, timer) \
47  if (verb) \
48  { \
49  std::cout << " - " << msg << ": " \
50  << timer.TimePerTest(1) << std::endl; \
51  }
52 
53 using namespace std;
54 using namespace Nektar::LibUtilities;
55 
56 namespace Nektar
57 {
58 namespace SpatialDomains
59 {
60 
61 /// Version of the Nektar++ HDF5 geometry format, which is embedded into the
62 /// main NEKTAR/GEOMETRY group as an attribute.
63 const unsigned int MeshGraphHDF5::FORMAT_VERSION = 1;
64 
65 std::string MeshGraphHDF5::className =
67  "HDF5", MeshGraphHDF5::create, "IO with HDF5 geometry");
68 
69 void MeshGraphHDF5::ReadGeometry(
70  DomainRangeShPtr rng,
71  bool fillGraph)
72 {
73  boost::ignore_unused(rng);
74 
75  ReadComposites();
76  ReadDomain();
77 
78  // Close up shop.
79  m_mesh->Close();
80  m_maps->Close();
81  m_file->Close();
82 
83  if (fillGraph)
84  {
85  MeshGraph::FillGraph();
86  }
87 }
88 
89 /**
90  * @brief Utility function to split a vector equally amongst a number of
91  * processors.
92  *
93  * @param vecsize Size of the total amount of work
94  * @param rank Rank of this process
95  * @param nprocs Number of processors in the group
96  *
97  * @return A pair with the offset this process should occupy, along with the
98  * count for the amount of work.
99  */
100 std::pair<size_t, size_t> SplitWork(size_t vecsize, int rank, int nprocs)
101 {
102  size_t div = vecsize / nprocs;
103  size_t rem = vecsize % nprocs;
104  if (rank < rem)
105  {
106  return std::make_pair(rank * (div + 1), div + 1);
107  }
108  else
109  {
110  return std::make_pair((rank - rem) * div + rem * (div + 1), div);
111  }
112 }
113 
114 template<class T, typename std::enable_if<T::kDim == 0, int>::type = 0>
115 inline int GetGeomDataDim(std::map<int, std::shared_ptr<T>> &geomMap)
116 {
117  boost::ignore_unused(geomMap);
118  return 3;
119 }
120 
121 template<class T, typename std::enable_if<T::kDim == 1, int>::type = 0>
122 inline int GetGeomDataDim(std::map<int, std::shared_ptr<T>> &geomMap)
123 {
124  boost::ignore_unused(geomMap);
125  return T::kNverts;
126 }
127 
128 template<class T, typename std::enable_if<T::kDim == 2, int>::type = 0>
129 inline int GetGeomDataDim(std::map<int, std::shared_ptr<T>> &geomMap)
130 {
131  boost::ignore_unused(geomMap);
132  return T::kNedges;
133 }
134 
135 template<class T, typename std::enable_if<T::kDim == 3, int>::type = 0>
136 inline int GetGeomDataDim(std::map<int, std::shared_ptr<T>> &geomMap)
137 {
138  boost::ignore_unused(geomMap);
139  return T::kNfaces;
140 }
141 
142 template<class ...T>
143 inline void UniqueValues(std::unordered_set<int> &unique)
144 {
145  boost::ignore_unused(unique);
146 }
147 
148 template<class ...T>
149 inline void UniqueValues(std::unordered_set<int> &unique,
150  const std::vector<int> &input,
151  T&... args)
152 {
153  for (auto i : input)
154  {
155  unique.insert(i);
156  }
157 
158  UniqueValues(unique, args...);
159 }
160 
161 /**
162  * @brief Partition the mesh
163  */
164 void MeshGraphHDF5::PartitionMesh(LibUtilities::SessionReaderSharedPtr session)
165 {
167  all.Start();
168  int err;
169  LibUtilities::CommSharedPtr comm = session->GetComm();
170  int rank = comm->GetRank(), nproc = comm->GetSize();
171 
172  // By default, only the root process will have read the session file, which
173  // is done to avoid every process needing to read the XML file. For HDF5, we
174  // don't care about this, so just have every process parse the session file.
175  if (rank > 0)
176  {
177  session->InitSession();
178  }
179 
180  // We use the XML geometry to find information about the HDF5 file.
181  m_session = session;
182  m_xmlGeom = m_session->GetElement("NEKTAR/GEOMETRY");
183  TiXmlAttribute *attr = m_xmlGeom->FirstAttribute();
184  m_meshPartitioned = true;
185  m_meshDimension = 3;
186  m_spaceDimension = 3;
187 
188  while (attr)
189  {
190  std::string attrName(attr->Name());
191  if (attrName == "DIM")
192  {
193  err = attr->QueryIntValue(&m_meshDimension);
194  ASSERTL0(err == TIXML_SUCCESS, "Unable to read mesh dimension.");
195  }
196  else if (attrName == "SPACE")
197  {
198  err = attr->QueryIntValue(&m_spaceDimension);
199  ASSERTL0(err == TIXML_SUCCESS, "Unable to read space dimension.");
200  }
201  else if (attrName == "PARTITION")
202  {
203  ASSERTL0(false,
204  "PARTITION parameter should only be used in XML meshes");
205  }
206  else if(attrName == "HDF5FILE")
207  {
208  m_hdf5Name = attr->Value();
209  ASSERTL1(err == TIXML_SUCCESS, "Unable to read hdf5 name.");
210  }
211  else if(attrName == "PARTITIONED")
212  {
213  ASSERTL0(false,
214  "PARTITIONED parameter should only be used in XML meshes");
215  }
216  else
217  {
218  std::string errstr("Unknown attribute: ");
219  errstr += attrName;
220  ASSERTL1(false, errstr.c_str());
221  }
222  // Get the next attribute.
223  attr = attr->Next();
224  }
225 
226  ASSERTL0(m_hdf5Name.size() > 0, "Unable to obtain mesh file name.");
227  ASSERTL0(m_meshDimension <= m_spaceDimension,
228  "Mesh dimension greater than space dimension.");
229 
230  // Open handle to the HDF5 mesh
231  LibUtilities::H5::PListSharedPtr parallelProps = H5::PList::Default();
232  m_readPL = H5::PList::Default();
233 
234  if (nproc > 1)
235  {
236  // Use MPI/O to access the file
237  parallelProps = H5::PList::FileAccess();
238  parallelProps->SetMpio(comm);
239  // Use collective IO
240  m_readPL = H5::PList::DatasetXfer();
241  m_readPL->SetDxMpioCollective();
242  }
243 
244  m_file = H5::File::Open(m_hdf5Name, H5F_ACC_RDONLY, parallelProps);
245 
246  auto root = m_file->OpenGroup("NEKTAR");
247  ASSERTL0(root, "Cannot find NEKTAR group in HDF5 file.");
248 
249  auto root2 = root->OpenGroup("GEOMETRY");
250  ASSERTL0(root2, "Cannot find NEKTAR/GEOMETRY group in HDF5 file.");
251 
252  // Check format version
253  unsigned int formatVersion;
254  H5::Group::AttrIterator attrIt = root2->attr_begin();
255  H5::Group::AttrIterator attrEnd = root2->attr_end();
256  for (; attrIt != attrEnd; ++attrIt)
257  {
258  if (*attrIt == "FORMAT_VERSION")
259  {
260  break;
261  }
262  }
263  ASSERTL0(attrIt != attrEnd,
264  "Unable to determine Nektar++ geometry HDF5 file version.");
265  root2->GetAttribute("FORMAT_VERSION", formatVersion);
266 
267  ASSERTL0(formatVersion <= FORMAT_VERSION,
268  "File format in " + m_hdf5Name + " is higher than supported in "
269  "this version of Nektar++");
270 
271  m_mesh = root2->OpenGroup("MESH");
272  ASSERTL0(m_mesh, "Cannot find NEKTAR/GEOMETRY/MESH group in HDF5 file.");
273  m_maps = root2->OpenGroup("MAPS");
274  ASSERTL0(m_mesh, "Cannot find NEKTAR/GEOMETRY/MAPS group in HDF5 file.");
275 
276  // Depending on dimension, read element IDs.
277  std::map<int, std::vector<std::tuple<
278  std::string, int, LibUtilities::ShapeType>>> dataSets;
279 
280  dataSets[1] = { make_tuple("SEG", 2, LibUtilities::eSegment) };
281  dataSets[2] = { make_tuple("TRI", 3, LibUtilities::eTriangle),
282  make_tuple("QUAD", 4, LibUtilities::eQuadrilateral) };
283  dataSets[3] = { make_tuple("TET", 4, LibUtilities::eTetrahedron),
284  make_tuple("PYR", 5, LibUtilities::ePyramid),
285  make_tuple("PRISM", 5, LibUtilities::ePrism),
286  make_tuple("HEX", 6, LibUtilities::eHexahedron) };
287 
288  // Read IDs for partitioning purposes
289  std::vector<MeshEntity> elmts;
290  std::vector<int> ids;
291 
292  // Map from element ID to 'row' which is a contiguous ordering required for
293  // parallel partitioning.
294  std::unordered_map<int, int> row2id, id2row;
295 
297  t.Start();
298  int rowCount = 0;
299  for (auto &it : dataSets[m_meshDimension])
300  {
301  std::string ds = std::get<0>(it);
302 
303  if (!m_mesh->ContainsDataSet(ds))
304  {
305  continue;
306  }
307 
308  // Open metadata dataset
309  H5::DataSetSharedPtr data = m_mesh->OpenDataSet(ds);
310  H5::DataSpaceSharedPtr space = data->GetSpace();
311  vector<hsize_t> dims = space->GetDims();
312 
313  H5::DataSetSharedPtr mdata = m_maps->OpenDataSet(ds);
314  H5::DataSpaceSharedPtr mspace = mdata->GetSpace();
315  vector<hsize_t> mdims = mspace->GetDims();
316 
317  // TODO: This could be done more intelligently; reads all IDs so that we
318  // can construct the dual graph of the mesh.
319  vector<int> tmpElmts, tmpIds;
320  mdata->Read(tmpIds, mspace, m_readPL);
321  data->Read(tmpElmts, space, m_readPL);
322 
323  const int nGeomData = std::get<1>(it);
324 
325  for (int i = 0, cnt = 0; i < tmpIds.size(); ++i, ++rowCount)
326  {
327  MeshEntity e;
328  row2id[rowCount] = tmpIds[i];
329  id2row[tmpIds[i]] = row2id[rowCount];
330  e.id = rowCount;
331  e.origId = tmpIds[i];
332  e.list = std::vector<unsigned int>(
333  &tmpElmts[cnt], &tmpElmts[cnt+nGeomData]);
334  elmts.push_back(e);
335  cnt += nGeomData;
336  }
337  }
338 
339  comm->Block();
340  t.Stop();
341 
342  bool verbRoot = rank == 0 &&
343  m_session->DefinesCmdLineArgument("verbose");
344 
345  if (verbRoot)
346  {
347  std::cout << "Reading HDF5 geometry..." << std::endl;
348  }
349 
350  TIME_RESULT(verbRoot, "initial read", t);
351 
352  // Check to see we have at least as many processors as elements.
353  size_t numElmt = elmts.size();
354  ASSERTL0(nproc <= numElmt,
355  "This mesh has more processors than elements!");
356 
357  // Calculate reasonably even distribution of processors for calling ptScotch
358  auto elRange = SplitWork(numElmt, rank, nproc);
359 
360  t.Start();
361 
362  // Construct map of element entities for partitioner.
363  std::map<int, MeshEntity> partElmts;
364  std::unordered_set<int> facetIDs;
365 
366  int vcnt = 0;
367 
368  for (int el = elRange.first; el < elRange.first + elRange.second;
369  ++el, ++vcnt)
370  {
371  MeshEntity elmt = elmts[el];
372  elmt.ghost = false;
373  partElmts[el] = elmt;
374 
375  for (auto &facet : elmt.list)
376  {
377  facetIDs.insert(facet);
378  }
379  }
380 
381  // Now identify ghost vertices for the graph. This could probably be
382  // improved.
383  int nLocal = vcnt;
384  for (int i = 0; i < numElmt; ++i)
385  {
386  // Ignore anything we already read.
387  if (i >= elRange.first && i < elRange.first + elRange.second)
388  {
389  //i += elRange.second - elRange.first - 1;
390  continue;
391  }
392 
393  MeshEntity elmt = elmts[i];
394  bool insert = false;
395 
396  // Check for connections to local elements.
397  for (auto &eId : elmt.list)
398  {
399  if (facetIDs.find(eId) != facetIDs.end())
400  {
401  insert = true;
402  break;
403  }
404  }
405 
406  if (insert)
407  {
408  elmt.ghost = true;
409  partElmts[elmt.id] = elmt;
410  }
411  }
412 
413  // Create partitioner. Default partitioner to use is PtScotch. Use ParMetis
414  // as default if it is installed. Override default with command-line flags
415  // if they are set.
416  string partitionerName = nproc > 1 ? "PtScotch" : "Scotch";
417  if (GetMeshPartitionFactory().ModuleExists("ParMetis"))
418  {
419  partitionerName = "ParMetis";
420  }
421  if (session->DefinesCmdLineArgument("use-parmetis"))
422  {
423  partitionerName = "ParMetis";
424  }
425  if (session->DefinesCmdLineArgument("use-ptscotch"))
426  {
427  partitionerName = "PtScotch";
428  }
429 
430  MeshPartitionSharedPtr partitioner =
432  partitionerName, session, m_meshDimension,
433  partElmts, CreateCompositeDescriptor(id2row));
434 
435  partitioner->PartitionMesh(nproc, true, false, nLocal);
436 
437  // Now we need to distribute vertex IDs to all the different processors.
438  t.Stop();
439  TIME_RESULT(verbRoot, "partitioning", t);
440 
441  // Each process now knows which rows of the dataset it needs to read for the
442  // elements of dimension m_meshDimension.
443  std::vector<unsigned int> tmp;
444  std::unordered_set<int> toRead;
445  partitioner->GetElementIDs(comm->GetRank(), tmp);
446 
447  for (auto &tmpId : tmp)
448  {
449  toRead.insert(row2id[tmpId]);
450  }
451 
452  // Since objects are going to be constructed starting from vertices, we now
453  // need to recurse down the geometry facet dimensions to figure out which
454  // rows to read from each dataset.
455  std::vector<int> vertIDs, segIDs, triIDs, quadIDs;
456  std::vector<int> tetIDs, prismIDs, pyrIDs, hexIDs;
457  std::vector<int> segData, triData, quadData, tetData;
458  std::vector<int> prismData, pyrData, hexData;
459  std::vector<NekDouble> vertData;
460 
461  if (m_meshDimension == 3)
462  {
463  t.Start();
464  // Read 3D data
465  ReadGeometryData(m_hexGeoms, "HEX", toRead, hexIDs, hexData);
466  ReadGeometryData(m_pyrGeoms, "PYR", toRead, pyrIDs, pyrData);
467  ReadGeometryData(m_prismGeoms, "PRISM", toRead, prismIDs, prismData);
468  ReadGeometryData(m_tetGeoms, "TET", toRead, tetIDs, tetData);
469 
470  toRead.clear();
471  UniqueValues(toRead, hexData, pyrData, prismData, tetData);
472  t.Stop();
473  TIME_RESULT(verbRoot, "read 3D elements", t);
474  }
475 
476  if (m_meshDimension >= 2)
477  {
478  t.Start();
479  // Read 2D data
480  ReadGeometryData(m_triGeoms, "TRI", toRead, triIDs, triData);
481  ReadGeometryData(m_quadGeoms, "QUAD", toRead, quadIDs, quadData);
482 
483  toRead.clear();
484  UniqueValues(toRead, triData, quadData);
485  t.Stop();
486  TIME_RESULT(verbRoot, "read 2D elements", t);
487  }
488 
489  if (m_meshDimension >= 1)
490  {
491  t.Start();
492  // Read 2D data
493  ReadGeometryData(m_segGeoms, "SEG", toRead, segIDs, segData);
494 
495  toRead.clear();
496  UniqueValues(toRead, segData);
497  t.Stop();
498  TIME_RESULT(verbRoot, "read 1D elements", t);
499  }
500 
501  t.Start();
502  ReadGeometryData(m_vertSet, "VERT", toRead, vertIDs, vertData);
503  t.Stop();
504  TIME_RESULT(verbRoot, "read 0D elements", t);
505 
506  // Now start to construct geometry objects, starting from vertices upwards.
507  t.Start();
508  FillGeomMap(m_vertSet, CurveMap(), vertIDs, vertData);
509  t.Stop();
510  TIME_RESULT(verbRoot, "construct 0D elements", t);
511 
512  if (m_meshDimension >= 1)
513  {
514  // Read curves
515  toRead.clear();
516  for (auto &edge : segIDs)
517  {
518  toRead.insert(edge);
519  }
520  ReadCurveMap(m_curvedEdges, "CURVE_EDGE", toRead);
521 
522  t.Start();
523  FillGeomMap(m_segGeoms, m_curvedEdges, segIDs, segData);
524  t.Stop();
525  TIME_RESULT(verbRoot, "construct 1D elements", t);
526  }
527 
528  if (m_meshDimension >= 2)
529  {
530  // Read curves
531  toRead.clear();
532  for (auto &face : triIDs)
533  {
534  toRead.insert(face);
535  }
536  for (auto &face : quadIDs)
537  {
538  toRead.insert(face);
539  }
540  ReadCurveMap(m_curvedFaces, "CURVE_FACE", toRead);
541 
542  t.Start();
543  FillGeomMap(m_triGeoms, m_curvedFaces, triIDs, triData);
544  FillGeomMap(m_quadGeoms, m_curvedFaces, quadIDs, quadData);
545  t.Stop();
546  TIME_RESULT(verbRoot, "construct 2D elements", t);
547  }
548 
549  if (m_meshDimension >= 3)
550  {
551  t.Start();
552  FillGeomMap(m_hexGeoms, CurveMap(), hexIDs, hexData);
553  FillGeomMap(m_prismGeoms, CurveMap(), prismIDs, prismData);
554  FillGeomMap(m_pyrGeoms, CurveMap(), pyrIDs, pyrData);
555  FillGeomMap(m_tetGeoms, CurveMap(), tetIDs, tetData);
556  t.Stop();
557  TIME_RESULT(verbRoot, "construct 3D elements", t);
558  }
559 
560  // Populate m_bndRegOrder.
561  if (m_session->DefinesElement("NEKTAR/CONDITIONS"))
562  {
563  std::set<int> vBndRegionIdList;
564  TiXmlElement *vConditions =
565  new TiXmlElement(*m_session->GetElement("Nektar/Conditions"));
566  TiXmlElement *vBndRegions =
567  vConditions->FirstChildElement("BOUNDARYREGIONS");
568  TiXmlElement *vItem;
569 
570  if (vBndRegions)
571  {
572  vItem = vBndRegions->FirstChildElement();
573  while (vItem)
574  {
575  std::string vSeqStr =
576  vItem->FirstChild()->ToText()->Value();
577  std::string::size_type indxBeg =
578  vSeqStr.find_first_of('[') + 1;
579  std::string::size_type indxEnd =
580  vSeqStr.find_last_of(']') - 1;
581  vSeqStr = vSeqStr.substr(indxBeg, indxEnd - indxBeg + 1);
582 
583  std::vector<unsigned int> vSeq;
584  ParseUtils::GenerateSeqVector(vSeqStr.c_str(), vSeq);
585 
586  int p = atoi(vItem->Attribute("ID"));
587  m_bndRegOrder[p] = vSeq;
588  vItem = vItem->NextSiblingElement();
589  }
590  }
591  }
592 
593  all.Stop();
594  TIME_RESULT(verbRoot, "total time", all);
595 }
596 
597 template<class T, typename DataType> void MeshGraphHDF5::ConstructGeomObject(
598  std::map<int, std::shared_ptr<T>> &geomMap, int id,
599  DataType *data, CurveSharedPtr curve)
600 {
601  boost::ignore_unused(geomMap, id, data, curve);
602 }
603 
604 template<> void MeshGraphHDF5::ConstructGeomObject(
605  std::map<int, std::shared_ptr<PointGeom>> &geomMap, int id,
606  NekDouble *data, CurveSharedPtr curve)
607 {
608  boost::ignore_unused(curve);
610  m_spaceDimension, id, data[0], data[1], data[2]);
611 }
612 
613 template<> void MeshGraphHDF5::ConstructGeomObject(
614  std::map<int, std::shared_ptr<SegGeom>> &geomMap, int id, int *data,
615  CurveSharedPtr curve)
616 {
617  PointGeomSharedPtr pts[2] = { GetVertex(data[0]), GetVertex(data[1]) };
619  id, m_spaceDimension, pts, curve);
620 }
621 
622 template<> void MeshGraphHDF5::ConstructGeomObject(
623  std::map<int, std::shared_ptr<TriGeom>> &geomMap, int id, int *data,
624  CurveSharedPtr curve)
625 {
626  SegGeomSharedPtr segs[3] = {
627  GetSegGeom(data[0]), GetSegGeom(data[1]), GetSegGeom(data[2]) };
628  geomMap[id] = MemoryManager<TriGeom>::AllocateSharedPtr(id, segs, curve);
629 }
630 
631 template<> void MeshGraphHDF5::ConstructGeomObject(
632  std::map<int, std::shared_ptr<QuadGeom>> &geomMap, int id, int *data,
633  CurveSharedPtr curve)
634 {
635  SegGeomSharedPtr segs[4] = {
636  GetSegGeom(data[0]), GetSegGeom(data[1]), GetSegGeom(data[2]),
637  GetSegGeom(data[3])
638  };
639  geomMap[id] = MemoryManager<QuadGeom>::AllocateSharedPtr(id, segs, curve);
640 }
641 
642 template<> void MeshGraphHDF5::ConstructGeomObject(
643  std::map<int, std::shared_ptr<TetGeom>> &geomMap, int id, int *data,
644  CurveSharedPtr curve)
645 {
646  boost::ignore_unused(curve);
647  TriGeomSharedPtr faces[4] = {
648  std::static_pointer_cast<TriGeom>(GetGeometry2D(data[0])),
649  std::static_pointer_cast<TriGeom>(GetGeometry2D(data[1])),
650  std::static_pointer_cast<TriGeom>(GetGeometry2D(data[2])),
651  std::static_pointer_cast<TriGeom>(GetGeometry2D(data[3]))
652  };
653 
654  auto tetGeom = MemoryManager<TetGeom>::AllocateSharedPtr(id, faces);
655  PopulateFaceToElMap(tetGeom, TetGeom::kNfaces);
656  geomMap[id] = tetGeom;
657 }
658 
659 template<> void MeshGraphHDF5::ConstructGeomObject(
660  std::map<int, std::shared_ptr<PyrGeom>> &geomMap, int id, int *data,
661  CurveSharedPtr curve)
662 {
663  boost::ignore_unused(curve);
664  Geometry2DSharedPtr faces[5] = {
665  GetGeometry2D(data[0]), GetGeometry2D(data[1]), GetGeometry2D(data[2]),
666  GetGeometry2D(data[3]), GetGeometry2D(data[4])
667  };
668 
669  auto pyrGeom = MemoryManager<PyrGeom>::AllocateSharedPtr(id, faces);
670  PopulateFaceToElMap(pyrGeom, PyrGeom::kNfaces);
671  geomMap[id] = pyrGeom;
672 }
673 
674 template<> void MeshGraphHDF5::ConstructGeomObject(
675  std::map<int, std::shared_ptr<PrismGeom>> &geomMap, int id, int *data,
676  CurveSharedPtr curve)
677 {
678  boost::ignore_unused(curve);
679  Geometry2DSharedPtr faces[5] = {
680  GetGeometry2D(data[0]), GetGeometry2D(data[1]), GetGeometry2D(data[2]),
681  GetGeometry2D(data[3]), GetGeometry2D(data[4])
682  };
683 
684  auto prismGeom = MemoryManager<PrismGeom>::AllocateSharedPtr(id, faces);
685  PopulateFaceToElMap(prismGeom, PrismGeom::kNfaces);
686  geomMap[id] = prismGeom;
687 }
688 
689 template<> void MeshGraphHDF5::ConstructGeomObject(
690  std::map<int, std::shared_ptr<HexGeom>> &geomMap, int id, int *data,
691  CurveSharedPtr curve)
692 {
693  boost::ignore_unused(curve);
694  QuadGeomSharedPtr faces[6] = {
695  std::static_pointer_cast<QuadGeom>(GetGeometry2D(data[0])),
696  std::static_pointer_cast<QuadGeom>(GetGeometry2D(data[1])),
697  std::static_pointer_cast<QuadGeom>(GetGeometry2D(data[2])),
698  std::static_pointer_cast<QuadGeom>(GetGeometry2D(data[3])),
699  std::static_pointer_cast<QuadGeom>(GetGeometry2D(data[4])),
700  std::static_pointer_cast<QuadGeom>(GetGeometry2D(data[5]))
701  };
702 
703  auto hexGeom = MemoryManager<HexGeom>::AllocateSharedPtr(id, faces);
704  PopulateFaceToElMap(hexGeom, HexGeom::kNfaces);
705  geomMap[id] = hexGeom;
706 }
707 
708 template<class T, typename DataType>
709 void MeshGraphHDF5::FillGeomMap(
710  std::map<int, std::shared_ptr<T>> &geomMap,
711  const CurveMap &curveMap,
712  std::vector<int> &ids,
713  std::vector<DataType> &geomData)
714 {
715  const int nGeomData = GetGeomDataDim(geomMap);
716  const int nRows = geomData.size() / nGeomData;
717  CurveSharedPtr empty;
718 
719  // Construct geometry object.
720  if (curveMap.size() > 0)
721  {
722  for(int i = 0, cnt = 0; i < nRows; i++, cnt += nGeomData)
723  {
724  auto cIt = curveMap.find(ids[i]);
725  ConstructGeomObject(
726  geomMap, ids[i], &geomData[cnt],
727  cIt == curveMap.end() ? empty : cIt->second);
728  }
729  }
730  else
731  {
732  for(int i = 0, cnt = 0; i < nRows; i++, cnt += nGeomData)
733  {
734  ConstructGeomObject(geomMap, ids[i], &geomData[cnt], empty);
735  }
736  }
737 }
738 
739 template<class T, typename DataType>
740 void MeshGraphHDF5::ReadGeometryData(
741  std::map<int, std::shared_ptr<T>> &geomMap,
742  std::string dataSet,
743  const std::unordered_set<int> &readIds,
744  std::vector<int> &ids,
745  std::vector<DataType> &geomData)
746 {
747  if (!m_mesh->ContainsDataSet(dataSet))
748  {
749  return;
750  }
751 
752  // Open mesh dataset
753  H5::DataSetSharedPtr data = m_mesh->OpenDataSet(dataSet);
754  H5::DataSpaceSharedPtr space = data->GetSpace();
755  vector<hsize_t> dims = space->GetDims();
756 
757  // Open metadata dataset
758  H5::DataSetSharedPtr mdata = m_maps->OpenDataSet(dataSet);
759  H5::DataSpaceSharedPtr mspace = mdata->GetSpace();
760  vector<hsize_t> mdims = mspace->GetDims();
761 
762  ASSERTL0(mdims[0] == dims[0], "map and data set lengths do not match");
763 
764  const int nGeomData = GetGeomDataDim(geomMap);
765 
766  // Read all IDs
767  vector<int> allIds;
768  mdata->Read(allIds, mspace);
769 
770  // Selective reading; clear data space range so that we can select certain
771  // rows from the datasets.
772  space->ClearRange();
773 
774  int i = 0;
775  std::vector<hsize_t> coords;
776  for (auto &id : allIds)
777  {
778  if (readIds.find(id) != readIds.end())
779  {
780  for (int j = 0; j < nGeomData; ++j)
781  {
782  coords.push_back(i);
783  coords.push_back(j);
784  }
785  ids.push_back(id);
786  }
787  ++i;
788  }
789 
790  space->SetSelection(coords.size() / 2, coords);
791 
792  // Read selected data.
793  data->Read(geomData, space, m_readPL);
794 }
795 
796 void MeshGraphHDF5::ReadCurveMap(
797  CurveMap &curveMap,
798  std::string dsName,
799  const std::unordered_set<int> &readIds)
800 {
801  // If dataset does not exist, exit.
802  if (!m_mesh->ContainsDataSet(dsName))
803  {
804  return;
805  }
806 
807  // Open up curve map data.
808  H5::DataSetSharedPtr curveData = m_mesh->OpenDataSet(dsName);
809  H5::DataSpaceSharedPtr curveSpace = curveData->GetSpace();
810 
811  // Open up ID data set.
812  H5::DataSetSharedPtr idData = m_maps->OpenDataSet(dsName);
813  H5::DataSpaceSharedPtr idSpace = idData->GetSpace();
814 
815  // Read all IDs and clear data space.
816  vector<int> ids, newIds;
817  idData->Read(ids, idSpace);
818  curveSpace->ClearRange();
819 
820  // Search IDs to figure out which curves to read.
821  vector<hsize_t> curveSel;
822 
823  int cnt = 0;
824  for (auto &id : ids)
825  {
826  if (readIds.find(id) != readIds.end())
827  {
828  curveSel.push_back(cnt);
829  curveSel.push_back(0);
830  curveSel.push_back(cnt);
831  curveSel.push_back(1);
832  curveSel.push_back(cnt);
833  curveSel.push_back(2);
834  newIds.push_back(id);
835  }
836 
837  ++cnt;
838  }
839 
840  // Check to see whether any processor will read anything
841  auto toRead = newIds.size();
842  m_session->GetComm()->AllReduce(toRead, LibUtilities::ReduceSum);
843 
844  if (toRead == 0)
845  {
846  return;
847  }
848 
849  // Now read curve map and read data.
850  vector<int> curveInfo;
851  curveSpace->SetSelection(curveSel.size() / 2, curveSel);
852  curveData->Read(curveInfo, curveSpace, m_readPL);
853 
854  curveSel.clear();
855 
856  std::unordered_map<int, int> curvePtOffset;
857 
858  // Construct curves. We'll populate nodes in a minute!
859  for (int i = 0, cnt = 0, cnt2 = 0; i < curveInfo.size() / 3; ++i, cnt += 3)
860  {
862  newIds[i], (LibUtilities::PointsType)curveInfo[cnt + 1]);
863 
864  curve->m_points.resize(curveInfo[cnt]);
865 
866  const int ptOffset = curveInfo[cnt + 2];
867 
868  for (int j = 0; j < curveInfo[cnt]; ++j)
869  {
870  // ptoffset gives us the row, multiply by 3 for number of
871  // coordinates.
872  curveSel.push_back(ptOffset + j);
873  curveSel.push_back(0);
874  curveSel.push_back(ptOffset + j);
875  curveSel.push_back(1);
876  curveSel.push_back(ptOffset + j);
877  curveSel.push_back(2);
878  }
879 
880  // Store the offset so we know to come back later on to fill in these
881  // points.
882  curvePtOffset[newIds[i]] = 3 * cnt2;
883  cnt2 += curveInfo[cnt];
884 
885  curveMap[newIds[i]] = curve;
886  }
887 
888  curveInfo.clear();
889 
890  // Open node data spacee.
891  H5::DataSetSharedPtr nodeData = m_mesh->OpenDataSet("CURVE_NODES");
892  H5::DataSpaceSharedPtr nodeSpace = nodeData->GetSpace();
893 
894  nodeSpace->ClearRange();
895  nodeSpace->SetSelection(curveSel.size() / 2, curveSel);
896 
897  vector<NekDouble> nodeRawData;
898  nodeData->Read(nodeRawData, nodeSpace, m_readPL);
899 
900  // Go back and populate data from nodes.
901  for (auto &cIt : curvePtOffset)
902  {
903  CurveSharedPtr curve = curveMap[cIt.first];
904 
905  // Create nodes.
906  int cnt = cIt.second;
907  for (int i = 0; i < curve->m_points.size(); ++i, cnt += 3)
908  {
909  curve->m_points[i] = MemoryManager<PointGeom>::AllocateSharedPtr(
910  0, m_spaceDimension, nodeRawData[cnt], nodeRawData[cnt+1],
911  nodeRawData[cnt+2]);
912  }
913  }
914 }
915 
916 void MeshGraphHDF5::ReadDomain()
917 {
918  map<int, CompositeSharedPtr> fullDomain;
919  H5::DataSetSharedPtr dst = m_mesh->OpenDataSet("DOMAIN");
920  H5::DataSpaceSharedPtr space = dst->GetSpace();
921 
922  vector<string> data;
923  dst->ReadVectorString(data, space, m_readPL);
924  GetCompositeList(data[0], fullDomain);
925  m_domain.push_back(fullDomain);
926 }
927 
928 void MeshGraphHDF5::ReadComposites()
929 {
930  string nm = "COMPOSITE";
931 
932  H5::DataSetSharedPtr data = m_mesh->OpenDataSet(nm);
933  H5::DataSpaceSharedPtr space = data->GetSpace();
934  vector<hsize_t> dims = space->GetDims();
935 
936  vector<string> comps;
937  data->ReadVectorString(comps, space);
938 
939  H5::DataSetSharedPtr mdata = m_maps->OpenDataSet(nm);
940  H5::DataSpaceSharedPtr mspace = mdata->GetSpace();
941  vector<hsize_t> mdims = mspace->GetDims();
942 
943  vector<int> ids;
944  mdata->Read(ids, mspace);
945 
946  for(int i = 0; i < dims[0]; i++)
947  {
948  string compStr = comps[i];
949 
950  char type;
951  istringstream strm(compStr);
952 
953  strm >> type;
954 
955  CompositeSharedPtr comp =
957 
958  string::size_type indxBeg = compStr.find_first_of('[') + 1;
959  string::size_type indxEnd = compStr.find_last_of(']') - 1;
960 
961  string indxStr = compStr.substr(indxBeg, indxEnd - indxBeg + 1);
962  vector<unsigned int> seqVector;
963 
964  ParseUtils::GenerateSeqVector(indxStr, seqVector);
965  m_compOrder[ids[i]] = seqVector;
966 
967  switch (type)
968  {
969  case 'V':
970  for (auto &i : seqVector)
971  {
972  auto it = m_vertSet.find(i);
973  if (it != m_vertSet.end())
974  {
975  comp->m_geomVec.push_back(it->second);
976  }
977  }
978  break;
979  case 'S':
980  case 'E':
981  for (auto &i : seqVector)
982  {
983  auto it = m_segGeoms.find(i);
984  if (it != m_segGeoms.end())
985  {
986  comp->m_geomVec.push_back(it->second);
987  }
988  }
989  break;
990  case 'Q':
991  for (auto &i : seqVector)
992  {
993  auto it = m_quadGeoms.find(i);
994  if (it != m_quadGeoms.end())
995  {
996  comp->m_geomVec.push_back(it->second);
997  }
998  }
999  break;
1000  case 'T':
1001  for (auto &i : seqVector)
1002  {
1003  auto it = m_triGeoms.find(i);
1004  if (it != m_triGeoms.end())
1005  {
1006  comp->m_geomVec.push_back(it->second);
1007  }
1008  }
1009  break;
1010  case 'F':
1011  for(auto &i : seqVector)
1012  {
1013  auto it1 = m_quadGeoms.find(i);
1014  if (it1 != m_quadGeoms.end())
1015  {
1016  comp->m_geomVec.push_back(it1->second);
1017  continue;
1018  }
1019  auto it2 = m_triGeoms.find(i);
1020  if (it2 != m_triGeoms.end())
1021  {
1022  comp->m_geomVec.push_back(it2->second);
1023  }
1024  }
1025  break;
1026  case 'A':
1027  for(auto &i : seqVector)
1028  {
1029  auto it = m_tetGeoms.find(i);
1030  if (it != m_tetGeoms.end())
1031  {
1032  comp->m_geomVec.push_back(it->second);
1033  }
1034  }
1035  break;
1036  case 'P':
1037  for(auto &i : seqVector)
1038  {
1039  auto it = m_pyrGeoms.find(i);
1040  if (it != m_pyrGeoms.end())
1041  {
1042  comp->m_geomVec.push_back(it->second);
1043  }
1044  }
1045  break;
1046  case 'R':
1047  for(auto &i : seqVector)
1048  {
1049  auto it = m_prismGeoms.find(i);
1050  if (it != m_prismGeoms.end())
1051  {
1052  comp->m_geomVec.push_back(it->second);
1053  }
1054  }
1055  break;
1056  case 'H':
1057  for(auto &i : seqVector)
1058  {
1059  auto it = m_hexGeoms.find(i);
1060  if (it != m_hexGeoms.end())
1061  {
1062  comp->m_geomVec.push_back(it->second);
1063  }
1064  }
1065  break;
1066  }
1067 
1068  if (comp->m_geomVec.size() > 0)
1069  {
1070  m_meshComposites[ids[i]] = comp;
1071  }
1072  }
1073 }
1074 
1075 CompositeDescriptor MeshGraphHDF5::CreateCompositeDescriptor(
1076  std::unordered_map<int, int> &id2row)
1077 {
1078  CompositeDescriptor ret;
1079 
1080  string nm = "COMPOSITE";
1081 
1082  H5::DataSetSharedPtr data = m_mesh->OpenDataSet(nm);
1083  H5::DataSpaceSharedPtr space = data->GetSpace();
1084  vector<hsize_t> dims = space->GetDims();
1085 
1086  vector<string> comps;
1087  data->ReadVectorString(comps, space);
1088 
1089  H5::DataSetSharedPtr mdata = m_maps->OpenDataSet(nm);
1090  H5::DataSpaceSharedPtr mspace = mdata->GetSpace();
1091  vector<hsize_t> mdims = mspace->GetDims();
1092 
1093  vector<int> ids;
1094  mdata->Read(ids, mspace);
1095 
1096  for (int i = 0; i < dims[0]; i++)
1097  {
1098  string compStr = comps[i];
1099 
1100  char type;
1101  istringstream strm(compStr);
1102 
1103  strm >> type;
1104 
1105  string::size_type indxBeg = compStr.find_first_of('[') + 1;
1106  string::size_type indxEnd = compStr.find_last_of(']') - 1;
1107 
1108  string indxStr = compStr.substr(indxBeg, indxEnd - indxBeg + 1);
1109  vector<unsigned int> seqVector;
1110  ParseUtils::GenerateSeqVector(indxStr, seqVector);
1111 
1112  LibUtilities::ShapeType shapeType;
1113 
1114  switch (type)
1115  {
1116  case 'V':
1117  shapeType = LibUtilities::ePoint;
1118  break;
1119  case 'S':
1120  case 'E':
1121  shapeType = LibUtilities::eSegment;
1122  break;
1123  case 'Q':
1124  case 'F':
1125  // Note that for HDF5, the composite descriptor is only used for
1126  // partitioning purposes so 'F' tag is not really going to be
1127  // critical in this context.
1128  shapeType = LibUtilities::eQuadrilateral;
1129  break;
1130  case 'T':
1131  shapeType = LibUtilities::eTriangle;
1132  break;
1133  case 'A':
1134  shapeType = LibUtilities::eTetrahedron;
1135  break;
1136  case 'P':
1137  shapeType = LibUtilities::ePyramid;
1138  break;
1139  case 'R':
1140  shapeType = LibUtilities::ePrism;
1141  break;
1142  case 'H':
1143  shapeType = LibUtilities::eHexahedron;
1144  break;
1145  }
1146 
1147  std::vector<int> filteredVector;
1148  for (auto &compElmt : seqVector)
1149  {
1150  if (id2row.find(compElmt) == id2row.end())
1151  {
1152  continue;
1153  }
1154 
1155  filteredVector.push_back(compElmt);
1156  }
1157 
1158  if (filteredVector.size() == 0)
1159  {
1160  continue;
1161  }
1162 
1163  ret[ids[i]] = std::make_pair(shapeType, filteredVector);
1164  }
1165 
1166  return ret;
1167 }
1168 
1169 
1170 template<class T, typename std::enable_if<T::kDim == 0, int>::type = 0>
1171 inline NekDouble GetGeomData(std::shared_ptr<T> &geom, int i)
1172 {
1173  return (*geom)(i);
1174 }
1175 
1176 template<class T, typename std::enable_if<T::kDim == 1, int>::type = 0>
1177 inline int GetGeomData(std::shared_ptr<T> &geom, int i)
1178 {
1179  return geom->GetVid(i);
1180 }
1181 
1182 template<class T, typename std::enable_if<T::kDim == 2, int>::type = 0>
1183 inline int GetGeomData(std::shared_ptr<T> &geom, int i)
1184 {
1185  return geom->GetEid(i);
1186 }
1187 
1188 template<class T, typename std::enable_if<T::kDim == 3, int>::type = 0>
1189 inline int GetGeomData(std::shared_ptr<T> &geom, int i)
1190 {
1191  return geom->GetFid(i);
1192 }
1193 
1194 template<class T>
1195 void MeshGraphHDF5::WriteGeometryMap(std::map<int, std::shared_ptr<T>> &geomMap,
1196  std::string datasetName)
1197 {
1198  typedef typename std::conditional<
1199  std::is_same<T, PointGeom>::value, NekDouble, int>::type DataType;
1200 
1201  const int nGeomData = GetGeomDataDim(geomMap);
1202  const size_t nGeom = geomMap.size();
1203 
1204  if (nGeom == 0)
1205  {
1206  return;
1207  }
1208 
1209  // Construct a map storing IDs
1210  vector<int> idMap(nGeom);
1211  vector<DataType> data(nGeom * nGeomData);
1212 
1213  int cnt1 = 0, cnt2 = 0;
1214  for (auto &it : geomMap)
1215  {
1216  idMap[cnt1++] = it.first;
1217 
1218  for (int j = 0; j < nGeomData; ++j)
1219  {
1220  data[cnt2 + j] = GetGeomData(it.second, j);
1221  }
1222 
1223  cnt2 += nGeomData;
1224  }
1225 
1226  vector<hsize_t> dims = { static_cast<hsize_t>(nGeom),
1227  static_cast<hsize_t>(nGeomData) };
1228  H5::DataTypeSharedPtr tp = H5::DataType::OfObject(data[0]);
1230  std::shared_ptr<H5::DataSpace>(new H5::DataSpace(dims));
1231  H5::DataSetSharedPtr dst = m_mesh->CreateDataSet(datasetName, tp, ds);
1232  dst->Write(data, ds);
1233 
1234  tp = H5::DataType::OfObject(idMap[0]);
1235  dims = { nGeom };
1236  ds = std::shared_ptr<H5::DataSpace>(new H5::DataSpace(dims));
1237  dst = m_maps->CreateDataSet(datasetName, tp, ds);
1238  dst->Write(idMap, ds);
1239 }
1240 
1241 void MeshGraphHDF5::WriteCurveMap(CurveMap &curves,
1242  std::string dsName,
1243  MeshCurvedPts &curvedPts,
1244  int &ptOffset,
1245  int &newIdx)
1246 {
1247  vector<int> data, map;
1248 
1249  // Compile curve data.
1250  for (auto &c : curves)
1251  {
1252  map.push_back(c.first);
1253  data.push_back(c.second->m_points.size());
1254  data.push_back(c.second->m_ptype);
1255  data.push_back(ptOffset);
1256 
1257  ptOffset += c.second->m_points.size();
1258 
1259  for (auto &pt : c.second->m_points)
1260  {
1261  MeshVertex v;
1262  v.id = newIdx;
1263  pt->GetCoords(v.x, v.y, v.z);
1264  curvedPts.pts.push_back(v);
1265  curvedPts.index.push_back(newIdx++);
1266  }
1267  }
1268 
1269  // Write data.
1270  vector<hsize_t> dims = { data.size() / 3, 3 };
1271  H5::DataTypeSharedPtr tp = H5::DataType::OfObject(data[0]);
1273  std::shared_ptr<H5::DataSpace>(new H5::DataSpace(dims));
1274  H5::DataSetSharedPtr dst = m_mesh->CreateDataSet(dsName, tp, ds);
1275  dst->Write(data, ds);
1276 
1277  tp = H5::DataType::OfObject(map[0]);
1278  dims = { map.size() };
1279  ds = std::shared_ptr<H5::DataSpace>(new H5::DataSpace(dims));
1280  dst = m_maps->CreateDataSet(dsName, tp, ds);
1281  dst->Write(map, ds);
1282 }
1283 
1284 void MeshGraphHDF5::WriteCurvePoints(MeshCurvedPts &curvedPts)
1285 {
1286  vector<double> vertData(curvedPts.pts.size() * 3);
1287 
1288  int cnt = 0;
1289  for (auto &pt : curvedPts.pts)
1290  {
1291  vertData[cnt++] = pt.x;
1292  vertData[cnt++] = pt.y;
1293  vertData[cnt++] = pt.z;
1294  }
1295 
1296  vector<hsize_t> dims = { curvedPts.pts.size(), 3 };
1297  H5::DataTypeSharedPtr tp = H5::DataType::OfObject(vertData[0]);
1299  std::shared_ptr<H5::DataSpace>(new H5::DataSpace(dims));
1300  H5::DataSetSharedPtr dst = m_mesh->CreateDataSet("CURVE_NODES", tp, ds);
1301  dst->Write(vertData, ds);
1302 }
1303 
1304 void MeshGraphHDF5::WriteComposites(CompositeMap &composites)
1305 {
1306  vector<string> comps;
1307 
1308  //dont need location map only a id map
1309  //will filter the composites per parition on read, its easier
1310  //composites do not need to be written in paralell.
1311  vector<int> c_map;
1312 
1313  for (auto &cIt : composites)
1314  {
1315  if (cIt.second->m_geomVec.size() == 0)
1316  {
1317  continue;
1318  }
1319 
1320  comps.push_back(GetCompositeString(cIt.second));
1321  c_map.push_back(cIt.first);
1322  }
1323 
1324  H5::DataTypeSharedPtr tp = H5::DataType::String();
1325  H5::DataSpaceSharedPtr ds = H5::DataSpace::OneD(comps.size());
1326  H5::DataSetSharedPtr dst = m_mesh->CreateDataSet("COMPOSITE", tp, ds);
1327  dst->WriteVectorString(comps, ds, tp);
1328 
1329  tp = H5::DataType::OfObject(c_map[0]);
1330  ds = H5::DataSpace::OneD(c_map.size());
1331  dst = m_maps->CreateDataSet("COMPOSITE", tp, ds);
1332  dst->Write(c_map, ds);
1333 }
1334 
1335 void MeshGraphHDF5::WriteDomain(vector<CompositeMap> &domain)
1336 {
1337  vector<unsigned int> idxList;
1338  for (auto cIt = domain[0].begin(); cIt != domain[0].end(); ++cIt)
1339  {
1340  idxList.push_back(cIt->first);
1341  }
1342  stringstream domString;
1343  vector<string> doms;
1344  doms.push_back(ParseUtils::GenerateSeqString(idxList));
1345 
1346  H5::DataTypeSharedPtr tp = H5::DataType::String();
1347  H5::DataSpaceSharedPtr ds = H5::DataSpace::OneD(doms.size());
1348  H5::DataSetSharedPtr dst = m_mesh->CreateDataSet("DOMAIN", tp, ds);
1349  dst->WriteVectorString(doms, ds, tp);
1350 }
1351 
1352 void MeshGraphHDF5::WriteGeometry(
1353  std::string &outfilename,
1354  bool defaultExp,
1355  const LibUtilities::FieldMetaDataMap &metadata)
1356 {
1357  boost::ignore_unused(metadata);
1358 
1359  vector<string> tmp;
1360  boost::split(tmp, outfilename, boost::is_any_of("."));
1361  string filenameXml = tmp[0] + ".xml";
1362  string filenameHdf5 = tmp[0] + ".nekg";
1363 
1364  //////////////////
1365  // XML part
1366  //////////////////
1367 
1368  //Check to see if a xml of the same name exists
1369  //if might have boundary conditions etc, we will just alter the geometry
1370  //tag if needed
1371  TiXmlDocument *doc = new TiXmlDocument;
1372  TiXmlElement *root;
1373  TiXmlElement *geomTag;
1374 
1375  if(boost::filesystem::exists(filenameXml.c_str()))
1376  {
1377  ifstream file(filenameXml.c_str());
1378  file >> (*doc);
1379  TiXmlHandle docHandle(doc);
1380  root = docHandle.FirstChildElement("NEKTAR").Element();
1381  ASSERTL0(root, "Unable to find NEKTAR tag in file.");
1382  geomTag = root->FirstChildElement("GEOMETRY");
1383  defaultExp = false;
1384  }
1385  else
1386  {
1387  TiXmlDeclaration *decl = new TiXmlDeclaration("1.0", "utf-8", "");
1388  doc->LinkEndChild(decl);
1389  root = new TiXmlElement("NEKTAR");
1390  doc->LinkEndChild(root);
1391 
1392  geomTag = new TiXmlElement("GEOMETRY");
1393  root->LinkEndChild(geomTag);
1394  }
1395 
1396  // Update attributes with dimensions.
1397  geomTag->SetAttribute("DIM", m_meshDimension);
1398  geomTag->SetAttribute("SPACE", m_spaceDimension);
1399  geomTag->SetAttribute("HDF5FILE", filenameHdf5);
1400 
1401  geomTag->Clear();
1402 
1403  if (defaultExp)
1404  {
1405  TiXmlElement *expTag = new TiXmlElement("EXPANSIONS");
1406 
1407  for (auto it = m_meshComposites.begin(); it != m_meshComposites.end();
1408  it++)
1409  {
1410  if (it->second->m_geomVec[0]->GetShapeDim() == m_meshDimension)
1411  {
1412  TiXmlElement *exp = new TiXmlElement("E");
1413  exp->SetAttribute(
1414  "COMPOSITE",
1415  "C[" + boost::lexical_cast<string>(it->first) + "]");
1416  exp->SetAttribute("NUMMODES", 4);
1417  exp->SetAttribute("TYPE", "MODIFIED");
1418  exp->SetAttribute("FIELDS", "u");
1419 
1420  expTag->LinkEndChild(exp);
1421  }
1422  }
1423  root->LinkEndChild(expTag);
1424  }
1425 
1426  doc->SaveFile(filenameXml);
1427 
1428  //////////////////
1429  // HDF5 part
1430  //////////////////
1431 
1432  // This is serial IO so we will just override any existing file.
1433  m_file = H5::File::Create(filenameHdf5, H5F_ACC_TRUNC);
1434  auto hdfRoot = m_file->CreateGroup("NEKTAR");
1435  auto hdfRoot2 = hdfRoot->CreateGroup("GEOMETRY");
1436 
1437  // Write format version.
1438  hdfRoot2->SetAttribute("FORMAT_VERSION", FORMAT_VERSION);
1439 
1440  // Create main groups.
1441  m_mesh = hdfRoot2->CreateGroup("MESH");
1442  m_maps = hdfRoot2->CreateGroup("MAPS");
1443 
1444  WriteGeometryMap(m_vertSet, "VERT");
1445  WriteGeometryMap(m_segGeoms, "SEG");
1446  if (m_meshDimension > 1)
1447  {
1448  WriteGeometryMap(m_triGeoms, "TRI");
1449  WriteGeometryMap(m_quadGeoms, "QUAD");
1450  }
1451  if (m_meshDimension > 2)
1452  {
1453  WriteGeometryMap(m_tetGeoms, "TET");
1454  WriteGeometryMap(m_pyrGeoms, "PYR");
1455  WriteGeometryMap(m_prismGeoms, "PRISM");
1456  WriteGeometryMap(m_hexGeoms, "HEX");
1457  }
1458 
1459  // Write curves
1460  int ptOffset = 0, newIdx = 0;
1461  MeshCurvedPts curvePts;
1462  WriteCurveMap(m_curvedEdges, "CURVE_EDGE", curvePts, ptOffset, newIdx);
1463  WriteCurveMap(m_curvedFaces, "CURVE_FACE", curvePts, ptOffset, newIdx);
1464  WriteCurvePoints(curvePts);
1465 
1466  // Write composites and domain.
1467  WriteComposites(m_meshComposites);
1468  WriteDomain(m_domain);
1469 }
1470 
1471 }
1472 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
#define TIME_RESULT(verb, msg, timer)
HDF5 DataSpace wrapper.
Definition: H5.h:315
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:199
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
std::shared_ptr< DataSpace > DataSpaceSharedPtr
Definition: H5.h:90
std::shared_ptr< PList > PListSharedPtr
Definition: H5.h:104
std::shared_ptr< DataType > DataTypeSharedPtr
Definition: H5.h:86
std::shared_ptr< DataSet > DataSetSharedPtr
Definition: H5.h:102
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:52
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:53
std::shared_ptr< QuadGeom > QuadGeomSharedPtr
Definition: HexGeom.h:46
void UniqueValues(std::unordered_set< int > &unique, const std::vector< int > &input, T &... args)
std::map< int, std::pair< LibUtilities::ShapeType, std::vector< int > > > CompositeDescriptor
Definition: MeshGraph.h:62
int GetGeomDataDim(std::map< int, std::shared_ptr< T >> &geomMap)
MeshPartitionFactory & GetMeshPartitionFactory()
std::shared_ptr< DomainRange > DomainRangeShPtr
Definition: MeshGraph.h:128
std::pair< size_t, size_t > SplitWork(size_t vecsize, int rank, int nprocs)
Utility function to split a vector equally amongst a number of processors.
std::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:61
std::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:62
std::map< int, CompositeSharedPtr > CompositeMap
Definition: MeshGraph.h:137
int GetGeomData(std::shared_ptr< T > &geom, int i)
std::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:62
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:59
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry.h:65
std::shared_ptr< MeshPartition > MeshPartitionSharedPtr
std::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:58
MeshGraphFactory & GetMeshGraphFactory()
Definition: MeshGraph.cpp:76
std::shared_ptr< Composite > CompositeSharedPtr
Definition: MeshGraph.h:136
double NekDouble
std::vector< MeshVertex > pts
mapping to access pts value.
std::vector< NekInt64 > index
id of this Point set
std::vector< unsigned int > list