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