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