Nektar++
Loading...
Searching...
No Matches
MeshGraphIOHDF5.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: MeshGraphIOHDF5.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 <boost/algorithm/string.hpp>
36#include <tinyxml.h>
37#include <type_traits>
38
45
46#define TIME_RESULT(verb, msg, timer) \
47 if (verb) \
48 { \
49 std::cout << " - " << msg << ": " << timer.TimePerTest(1) << "\n" \
50 << std::endl; \
51 }
52
53using namespace Nektar::LibUtilities;
54
56{
57
58/// Version of the Nektar++ HDF5 geometry format, which is embedded into the
59/// main NEKTAR/GEOMETRY group as an attribute.
60const unsigned int MeshGraphIOHDF5::FORMAT_VERSION = 2;
61
64 "HDF5", MeshGraphIOHDF5::create, "IO with HDF5 geometry");
65
67{
69 ReadDomain();
70
71 m_meshGraph->ReadExpansionInfo();
72
73 // Close up shop.
74 m_mesh->Close();
75 m_maps->Close();
76 m_file->Close();
77
78 if (fillGraph)
79 {
80 m_meshGraph->FillGraph();
81 }
82}
83
84/**
85 * @brief Utility function to split a vector equally amongst a number of
86 * processors.
87 *
88 * @param vecsize Size of the total amount of work
89 * @param rank Rank of this process
90 * @param nprocs Number of processors in the group
91 *
92 * @return A pair with the offset this process should occupy, along with the
93 * count for the amount of work.
94 */
95std::pair<size_t, size_t> SplitWork(size_t vecsize, int rank, int nprocs)
96{
97 size_t div = vecsize / nprocs;
98 size_t rem = vecsize % nprocs;
99 if (rank < rem)
100 {
101 return std::make_pair(rank * (div + 1), div + 1);
102 }
103 else
104 {
105 return std::make_pair((rank - rem) * div + rem * (div + 1), div);
106 }
107}
108
109template <class T, typename std::enable_if<T::kDim == 0, int>::type = 0>
110inline int GetGeomDataDim([[maybe_unused]] GeomMapView<T> &geomMap)
111{
112 return 3;
113}
114
115template <class T, typename std::enable_if<T::kDim == 1, int>::type = 0>
116inline int GetGeomDataDim([[maybe_unused]] GeomMapView<T> &geomMap)
117{
118 return T::kNverts;
119}
120
121template <class T, typename std::enable_if<T::kDim == 2, int>::type = 0>
122inline int GetGeomDataDim([[maybe_unused]] GeomMapView<T> &geomMap)
123{
124 return T::kNedges;
125}
126
127template <class T, typename std::enable_if<T::kDim == 3, int>::type = 0>
128inline int GetGeomDataDim([[maybe_unused]] GeomMapView<T> &geomMap)
129{
130 return T::kNfaces;
131}
132
133template <class... T>
134inline void UniqueValues([[maybe_unused]] std::unordered_set<int> &unique)
135{
136}
137
138template <class... T>
139inline void UniqueValues(std::unordered_set<int> &unique,
140 const std::vector<int> &input, T &...args)
141{
142 for (auto i : input)
143 {
144 unique.insert(i);
145 }
146
147 UniqueValues(unique, args...);
148}
149
150std::string MeshGraphIOHDF5::cmdSwitch =
152 "use-hdf5-node-comm", "",
153 "Use a per-node communicator for HDF5 partitioning.");
154
155/**
156 * @brief Partition the mesh
157 */
160
161{
163 all.Start();
164 int err;
165 LibUtilities::CommSharedPtr comm = session->GetComm();
166 LibUtilities::CommSharedPtr commMesh = comm->GetRowComm();
167 const bool isRoot = comm->TreatAsRankZero();
168
169 // By default, only the root process will have read the session file, which
170 // is done to avoid every process needing to read the XML file. For HDF5, we
171 // don't care about this, so just have every process parse the session file.
172 if (!isRoot)
173 {
174 session->InitSession();
175 }
176
177 // We use the XML geometry to find information about the HDF5 file.
178 m_session = session;
179 m_xmlGeom = session->GetElement("NEKTAR/GEOMETRY");
180 TiXmlAttribute *attr = m_xmlGeom->FirstAttribute();
181 int meshDimension = 3;
182 int spaceDimension = 3;
183
184 while (attr)
185 {
186 std::string attrName(attr->Name());
187 if (attrName == "DIM")
188 {
189 err = attr->QueryIntValue(&meshDimension);
190 ASSERTL0(err == TIXML_SUCCESS, "Unable to read mesh dimension.");
191 }
192 else if (attrName == "SPACE")
193 {
194 err = attr->QueryIntValue(&spaceDimension);
195 ASSERTL0(err == TIXML_SUCCESS, "Unable to read space dimension.");
196 }
197 else if (attrName == "PARTITION")
198 {
199 ASSERTL0(false,
200 "PARTITION parameter should only be used in XML meshes");
201 }
202 else if (attrName == "HDF5FILE")
203 {
204 m_hdf5Name = attr->Value();
205 }
206 else if (attrName == "PARTITIONED")
207 {
208 ASSERTL0(false,
209 "PARTITIONED parameter should only be used in XML meshes");
210 }
211 else
212 {
213 std::string errstr("Unknown attribute: ");
214 errstr += attrName;
215 ASSERTL1(false, errstr.c_str());
216 }
217 // Get the next attribute.
218 attr = attr->Next();
219 }
220
221 ASSERTL0(m_hdf5Name.size() > 0, "Unable to obtain mesh file name.");
222 ASSERTL0(meshDimension <= spaceDimension,
223 "Mesh dimension greater than space dimension.");
224
225 m_meshGraph->SetMeshDimension(meshDimension);
226 m_meshGraph->SetSpaceDimension(spaceDimension);
227
228 // Open handle to the HDF5 mesh
231
232 if (commMesh->GetSize() > 1)
233 {
234 // Use MPI/O to access the file
235 parallelProps = H5::PList::FileAccess();
236 parallelProps->SetMpio(commMesh);
237 // Use collective IO
239 m_readPL->SetDxMpioCollective();
240 }
241
242 m_file = H5::File::Open(m_hdf5Name, H5F_ACC_RDONLY, parallelProps);
243
244 auto root = m_file->OpenGroup("NEKTAR");
245 ASSERTL0(root, "Cannot find NEKTAR group in HDF5 file.");
246
247 auto root2 = root->OpenGroup("GEOMETRY");
248 ASSERTL0(root2, "Cannot find NEKTAR/GEOMETRY group in HDF5 file.");
249
250 // Check format version
251 H5::Group::AttrIterator attrIt = root2->attr_begin();
252 H5::Group::AttrIterator attrEnd = root2->attr_end();
253 for (; attrIt != attrEnd; ++attrIt)
254 {
255 if (*attrIt == "FORMAT_VERSION")
256 {
257 break;
258 }
259 }
260 ASSERTL0(attrIt != attrEnd,
261 "Unable to determine Nektar++ geometry HDF5 file version.");
262 root2->GetAttribute("FORMAT_VERSION", m_inFormatVersion);
263
265 "File format in " + m_hdf5Name +
266 " is higher than supported in "
267 "this version of Nektar++");
268
269 m_mesh = root2->OpenGroup("MESH");
270 ASSERTL0(m_mesh, "Cannot find NEKTAR/GEOMETRY/MESH group in HDF5 file.");
271 m_maps = root2->OpenGroup("MAPS");
272 ASSERTL0(m_mesh, "Cannot find NEKTAR/GEOMETRY/MAPS group in HDF5 file.");
273
275 {
276 return;
277 }
278
279 if (m_meshGraph->GetDomainRange() &&
280 m_meshGraph->GetDomainRange()->m_compElmts)
281 {
282 m_meshGraph->GetDomainRange()->m_compElmts = meshDimension;
283 }
284
285 SetupCompositeRange(m_meshGraph->GetDomainRange());
286
287 m_meshPartitioned = true;
288 m_meshGraph->SetMeshPartitioned(true);
289
290 // Depending on dimension, read element IDs.
291 std::map<int,
292 std::vector<std::tuple<std::string, int, LibUtilities::ShapeType>>>
293 dataSets;
294
295 dataSets[1] = {{"SEG", 2, LibUtilities::eSegment}};
296 dataSets[2] = {{"TRI", 3, LibUtilities::eTriangle},
297 {"QUAD", 4, LibUtilities::eQuadrilateral}};
298 dataSets[3] = {{"TET", 4, LibUtilities::eTetrahedron},
299 {"PYR", 5, LibUtilities::ePyramid},
300 {"PRISM", 5, LibUtilities::ePrism},
301 {"HEX", 6, LibUtilities::eHexahedron}};
302
303 const bool verbRoot = isRoot && session->DefinesCmdLineArgument("verbose");
304
305 if (verbRoot)
306 {
307 std::cout << "Reading HDF5 geometry..." << std::endl;
308 }
309
310 // If we want to use an multi-level communicator, then split the
311 // communicator at this point. We set by default the inter-node communicator
312 // to be the normal communicator: this way if the multi-level partitioning
313 // is disabled we proceed as 'normal'.
314 LibUtilities::CommSharedPtr innerComm, interComm = comm;
315 int innerRank = 0, innerSize = 1,
316 interRank = interComm->GetRowComm()->GetRank(),
317 interSize = interComm->GetRowComm()->GetSize();
318
319 if (session->DefinesCmdLineArgument("use-hdf5-node-comm"))
320 {
321 ASSERTL0(comm->GetSize() == commMesh->GetSize(),
322 "--use-hdf5-node-comm not available with Parallel-in-Time")
323
324 auto splitComm = comm->SplitCommNode();
325 innerComm = splitComm.first;
326 interComm = splitComm.second;
327 innerRank = innerComm->GetRank();
328 innerSize = innerComm->GetSize();
329
330 if (innerRank == 0)
331 {
332 interRank = interComm->GetRank();
333 interSize = interComm->GetSize();
334 }
335 }
336
337 // Unordered set of rows of the dataset this process needs to read for the
338 // elements of dimension meshDimension.
339 std::unordered_set<int> toRead;
340
341 // Calculate reasonably even distribution of processors for calling
342 // ptScotch. We'll do this work on only one process per node.
344 t.Start();
345
346 if (verbRoot)
347 {
348 std::cout << " - beginning partitioning" << std::endl;
349 }
350
351 // Perform initial read if either (a) we are on all ranks and multi-level
352 // partitioning is not enabled; or (b) rank 0 of all nodes if it is.
353 if (innerRank == 0)
354 {
356 t2.Start();
357
358 const bool verbRoot2 =
359 isRoot && session->DefinesCmdLineArgument("verbose");
360
361 // Read IDs for partitioning purposes
362 std::vector<int> ids;
363
364 // Map from element ID to 'row' which is a contiguous ordering required
365 // for parallel partitioning.
366 std::vector<MeshEntity> elmts;
367 std::unordered_map<int, int> row2id, id2row;
368
372
373 if (innerComm)
374 {
375 // For per-node partitioning, create a temporary reader (otherwise
376 // communicators are inconsistent).
377 auto parallelProps = H5::PList::FileAccess();
378 parallelProps->SetMpio(interComm);
379
380 // Use collective IO
381 readPL = H5::PList::DatasetXfer();
382 readPL->SetDxMpioCollective();
383 file = H5::File::Open(m_hdf5Name, H5F_ACC_RDONLY, parallelProps);
384
385 auto root = file->OpenGroup("NEKTAR");
386 auto root2 = root->OpenGroup("GEOMETRY");
387 mesh = root2->OpenGroup("MESH");
388 maps = root2->OpenGroup("MAPS");
389 }
390
391 int rowCount = 0;
392 for (auto &it : dataSets[meshDimension])
393 {
394 std::string ds = std::get<0>(it);
395
396 if (!mesh->ContainsDataSet(ds))
397 {
398 continue;
399 }
400
401 // Open metadata dataset
402 H5::DataSetSharedPtr data = mesh->OpenDataSet(ds);
403 H5::DataSpaceSharedPtr space = data->GetSpace();
404 std::vector<hsize_t> dims = space->GetDims();
405
406 H5::DataSetSharedPtr mdata = maps->OpenDataSet(ds);
407 H5::DataSpaceSharedPtr mspace = mdata->GetSpace();
408 std::vector<hsize_t> mdims = mspace->GetDims();
409
410 // TODO: This could perhaps be done more intelligently; reads all
411 // IDs for the top-level elements so that we can construct the dual
412 // graph of the mesh.
413 std::vector<int> tmpElmts, tmpIds;
414 mdata->Read(tmpIds, mspace, readPL);
415 data->Read(tmpElmts, space, readPL);
416
417 const int nGeomData = std::get<1>(it);
418
419 if (m_meshGraph->GetDomainRange() ==
421 {
422 // avoid range checking on larger meshes if not required
423 for (int i = 0, cnt = 0; i < tmpIds.size(); ++i, ++rowCount)
424 {
425 MeshEntity e;
426 row2id[rowCount] = tmpIds[i];
427 id2row[tmpIds[i]] = row2id[rowCount];
428 e.id = rowCount;
429 e.origId = tmpIds[i];
430 e.ghost = false;
431 e.list = std::vector<unsigned int>(
432 &tmpElmts[cnt], &tmpElmts[cnt + nGeomData]);
433 elmts.push_back(e);
434 cnt += nGeomData;
435 }
436 }
437 else
438 {
439 // avoid range checking on larger meshes if not required
440 for (int i = 0, cnt = 0; i < tmpIds.size(); ++i, ++rowCount)
441 {
442 MeshEntity e;
443 row2id[rowCount] = tmpIds[i];
444 id2row[tmpIds[i]] = row2id[rowCount];
445 e.id = rowCount;
446 e.origId = tmpIds[i];
447 e.ghost = false;
448 e.list = std::vector<unsigned int>(
449 &tmpElmts[cnt], &tmpElmts[cnt + nGeomData]);
450 if (m_meshGraph->CheckRange(e))
451 {
452 elmts.push_back(e);
453 }
454 cnt += nGeomData;
455 }
456 }
457 }
458
459 interComm->GetRowComm()->Block();
460
461 t2.Stop();
462 TIME_RESULT(verbRoot2, " - initial read", t2);
463 t2.Start();
464
465 // Do not partition in serial since postprocessing may not
466 // lead to very suitable element distribution that then leads
467 // to challenges with Scotch
468 if (commMesh->GetSize() > 1)
469 {
470 // Check to see we have at least as many processors as elements.
471 size_t numElmt = elmts.size();
472 ASSERTL0(commMesh->GetSize() <= numElmt,
473 "This mesh has more processors than elements!");
474
475 auto elRange = SplitWork(numElmt, interRank, interSize);
476
477 // Construct map of element entities for partitioner.
478 std::map<int, MeshEntity> partElmts;
479 std::unordered_set<int> facetIDs;
480
481 int vcnt = 0;
482
483 for (int el = elRange.first; el < elRange.first + elRange.second;
484 ++el, ++vcnt)
485 {
486 MeshEntity elmt = elmts[el];
487 elmt.ghost = false;
488 partElmts[elmt.id] = elmt;
489
490 for (auto &facet : elmt.list)
491 {
492 facetIDs.insert(facet);
493 }
494 }
495
496 // Now identify ghost vertices for the graph. This could also
497 // probably be improved.
498 int nLocal = vcnt;
499 for (int i = 0; i < numElmt; ++i)
500 {
501 // Ignore anything we already read.
502 if (i >= elRange.first && i < elRange.first + elRange.second)
503 {
504 continue;
505 }
506
507 MeshEntity elmt = elmts[i];
508 bool insert = false;
509
510 // Check for connections to local elements.
511 for (auto &eId : elmt.list)
512 {
513 if (facetIDs.find(eId) != facetIDs.end())
514 {
515 insert = true;
516 break;
517 }
518 }
519
520 if (insert)
521 {
522 elmt.ghost = true;
523 partElmts[elmt.id] = elmt;
524 }
525 }
526
527 // Create partitioner. Default partitioner to use is PtScotch. Use
528 // ParMetis as default if it is installed. Override default with
529 // command-line flags if they are set.
530 std::string partitionerName =
531 commMesh->GetSize() > 1 ? "PtScotch" : "Scotch";
532 if (GetMeshPartitionFactory().ModuleExists("ParMetis"))
533 {
534 partitionerName = "ParMetis";
535 }
536 if (session->DefinesCmdLineArgument("use-parmetis"))
537 {
538 partitionerName = "ParMetis";
539 }
540 if (session->DefinesCmdLineArgument("use-ptscotch"))
541 {
542 partitionerName = "PtScotch";
543 }
544
545 MeshPartitionSharedPtr partitioner =
547 partitionerName, session, interComm, meshDimension,
548 partElmts, CreateCompositeDescriptor(id2row));
549
550 t2.Stop();
551 TIME_RESULT(verbRoot2, " - partitioner setup", t2);
552 t2.Start();
553
554 partitioner->PartitionMesh(interSize, true, false, nLocal);
555 t2.Stop();
556 TIME_RESULT(verbRoot2, " - partitioning", t2);
557 t2.Start();
558
559 // Now construct a second graph that is partitioned in serial by
560 // this rank.
561 std::vector<unsigned int> nodeElmts;
562 partitioner->GetElementIDs(interRank, nodeElmts);
563
564 if (innerSize > 1)
565 {
566 // Construct map of element entities for partitioner.
567 std::map<int, MeshEntity> partElmts;
568 std::unordered_map<int, int> row2elmtid, elmtid2row;
569
570 int vcnt = 0;
571
572 // We need to keep track of which elements in the new partition
573 // correspond to elemental IDs for later (in a similar manner to
574 // row2id).
575 for (auto &elmtRow : nodeElmts)
576 {
577 row2elmtid[vcnt] = elmts[elmtRow].origId;
578 elmtid2row[elmts[elmtRow].origId] = vcnt;
579 MeshEntity elmt = elmts[elmtRow];
580 elmt.ghost = false;
581 partElmts[vcnt++] = elmt;
582 }
583
584 // Create temporary serial communicator for serial partitioning.
586 "Serial", 0, 0);
587
588 MeshPartitionSharedPtr partitioner =
590 "Scotch", session, tmpComm, meshDimension, partElmts,
591 CreateCompositeDescriptor(elmtid2row));
592
593 t2.Stop();
594 TIME_RESULT(verbRoot2, " - inner partition setup", t2);
595 t2.Start();
596
597 partitioner->PartitionMesh(innerSize, true, false, 0);
598
599 t2.Stop();
600 TIME_RESULT(verbRoot2, " - inner partitioning", t2);
601 t2.Start();
602
603 // Send contributions to remaining processors.
604 for (int i = 1; i < innerSize; ++i)
605 {
606 std::vector<unsigned int> tmp;
607 partitioner->GetElementIDs(i, tmp);
608 size_t tmpsize = tmp.size();
609 for (int j = 0; j < tmpsize; ++j)
610 {
611 tmp[j] = row2elmtid[tmp[j]];
612 }
613 innerComm->Send(i, tmpsize);
614 innerComm->Send(i, tmp);
615 }
616
617 t2.Stop();
618 TIME_RESULT(verbRoot2, " - inner partition scatter", t2);
619
620 std::vector<unsigned int> tmp;
621 partitioner->GetElementIDs(0, tmp);
622
623 for (auto &tmpId : tmp)
624 {
625 toRead.insert(row2elmtid[tmpId]);
626 }
627 }
628 else
629 {
630 for (auto &tmpId : nodeElmts)
631 {
632 toRead.insert(row2id[tmpId]);
633 }
634 }
635 }
636 else // Serial: Fill toRead with Elmt.origId
637 {
638 for (auto &tmpId : elmts)
639 {
640 toRead.insert(tmpId.origId);
641 }
642 }
643 }
644 else
645 {
646 // For multi-level partitioning, the innermost rank receives its
647 // partitions from rank 0 on each node.
648 size_t tmpSize;
649 innerComm->Recv(0, tmpSize);
650 std::vector<unsigned int> tmp(tmpSize);
651 innerComm->Recv(0, tmp);
652
653 for (auto &tmpId : tmp)
654 {
655 toRead.insert(tmpId);
656 }
657 }
658
659 t.Stop();
660 TIME_RESULT(verbRoot, "partitioning total", t);
661
662 // Since objects are going to be constructed starting from vertices, we
663 // now need to recurse down the geometry facet dimensions to figure out
664 // which rows to read from each dataset.
665 std::vector<int> vertIDs, segIDs, triIDs, quadIDs;
666 std::vector<int> tetIDs, prismIDs, pyrIDs, hexIDs;
667 std::vector<int> segData, triData, quadData, tetData;
668 std::vector<int> prismData, pyrData, hexData;
669 std::vector<NekDouble> vertData;
670
671 auto &vertSet = m_meshGraph->GetGeomMap<PointGeom>();
672 auto &segGeoms = m_meshGraph->GetGeomMap<SegGeom>();
673 auto &triGeoms = m_meshGraph->GetGeomMap<TriGeom>();
674 auto &quadGeoms = m_meshGraph->GetGeomMap<QuadGeom>();
675 auto &hexGeoms = m_meshGraph->GetGeomMap<HexGeom>();
676 auto &pyrGeoms = m_meshGraph->GetGeomMap<PyrGeom>();
677 auto &prismGeoms = m_meshGraph->GetGeomMap<PrismGeom>();
678 auto &tetGeoms = m_meshGraph->GetGeomMap<TetGeom>();
679 auto &curvedEdges = m_meshGraph->GetCurvedEdges();
680 auto &curvedFaces = m_meshGraph->GetCurvedFaces();
681
682 if (meshDimension == 3)
683 {
684 t.Start();
685 // Read 3D data
686 ReadGeometryData(hexGeoms, "HEX", toRead, hexIDs, hexData);
687 ReadGeometryData(pyrGeoms, "PYR", toRead, pyrIDs, pyrData);
688 ReadGeometryData(prismGeoms, "PRISM", toRead, prismIDs, prismData);
689 ReadGeometryData(tetGeoms, "TET", toRead, tetIDs, tetData);
690
691 toRead.clear();
692 UniqueValues(toRead, hexData, pyrData, prismData, tetData);
693 t.Stop();
694 TIME_RESULT(verbRoot, "read 3D elements", t);
695 }
696
697 if (meshDimension >= 2)
698 {
699 t.Start();
700 // Read 2D data
701 ReadGeometryData(triGeoms, "TRI", toRead, triIDs, triData);
702 ReadGeometryData(quadGeoms, "QUAD", toRead, quadIDs, quadData);
703
704 toRead.clear();
705 UniqueValues(toRead, triData, quadData);
706 t.Stop();
707 TIME_RESULT(verbRoot, "read 2D elements", t);
708 }
709
710 if (meshDimension >= 1)
711 {
712 t.Start();
713 // Read 1D data
714 ReadGeometryData(segGeoms, "SEG", toRead, segIDs, segData);
715
716 toRead.clear();
717 UniqueValues(toRead, segData);
718 t.Stop();
719 TIME_RESULT(verbRoot, "read 1D elements", t);
720 }
721
722 t.Start();
723 ReadGeometryData(vertSet, "VERT", toRead, vertIDs, vertData);
724 t.Stop();
725 TIME_RESULT(verbRoot, "read 0D elements", t);
726
727 // Now start to construct geometry objects, starting from vertices
728 // upwards.
729 t.Start();
730 FillGeomMap(vertSet, CurveMap(), vertIDs, vertData);
731 t.Stop();
732 TIME_RESULT(verbRoot, "construct 0D elements", t);
733
734 if (meshDimension >= 1)
735 {
736 // Read curves
737 toRead.clear();
738 for (auto &edge : segIDs)
739 {
740 toRead.insert(edge);
741 }
742 ReadCurveMap(curvedEdges, "CURVE_EDGE", toRead);
743
744 t.Start();
745 FillGeomMap(segGeoms, curvedEdges, segIDs, segData);
746 t.Stop();
747 TIME_RESULT(verbRoot, "construct 1D elements", t);
748 }
749
750 if (meshDimension >= 2)
751 {
752 // Read curves
753 toRead.clear();
754 for (auto &face : triIDs)
755 {
756 toRead.insert(face);
757 }
758 for (auto &face : quadIDs)
759 {
760 toRead.insert(face);
761 }
762 ReadCurveMap(curvedFaces, "CURVE_FACE", toRead);
763
764 t.Start();
765 FillGeomMap(triGeoms, curvedFaces, triIDs, triData);
766 FillGeomMap(quadGeoms, curvedFaces, quadIDs, quadData);
767 t.Stop();
768 TIME_RESULT(verbRoot, "construct 2D elements", t);
769 }
770
771 if (meshDimension >= 3)
772 {
773 t.Start();
774 FillGeomMap(hexGeoms, CurveMap(), hexIDs, hexData);
775 FillGeomMap(prismGeoms, CurveMap(), prismIDs, prismData);
776 FillGeomMap(pyrGeoms, CurveMap(), pyrIDs, pyrData);
777 FillGeomMap(tetGeoms, CurveMap(), tetIDs, tetData);
778 t.Stop();
779 TIME_RESULT(verbRoot, "construct 3D elements", t);
780 }
781
782 // Populate m_bndRegOrder.
783 if (session->DefinesElement("NEKTAR/CONDITIONS"))
784 {
785 std::set<int> vBndRegionIdList;
786 TiXmlElement *vConditions =
787 new TiXmlElement(*session->GetElement("Nektar/Conditions"));
788 TiXmlElement *vBndRegions =
789 vConditions->FirstChildElement("BOUNDARYREGIONS");
790 // Use fine-level for mesh partition (Parallel-in-Time)
792 TiXmlElement *vItem;
793
794 if (vBndRegions)
795 {
796 auto &graph_bndRegOrder = m_meshGraph->GetBndRegionOrdering();
797 vItem = vBndRegions->FirstChildElement();
798 while (vItem)
799 {
800 std::string vSeqStr = vItem->FirstChild()->ToText()->Value();
801 std::string::size_type indxBeg = vSeqStr.find_first_of('[') + 1;
802 std::string::size_type indxEnd = vSeqStr.find_last_of(']') - 1;
803 vSeqStr = vSeqStr.substr(indxBeg, indxEnd - indxBeg + 1);
804
805 std::vector<unsigned int> vSeq;
806 ParseUtils::GenerateSeqVector(vSeqStr.c_str(), vSeq);
807
808 int p = atoi(vItem->Attribute("ID"));
809 m_bndRegOrder[p] = vSeq;
810 graph_bndRegOrder[p] = vSeq;
811 vItem = vItem->NextSiblingElement();
812 }
813 }
814 }
815
816 all.Stop();
817 TIME_RESULT(verbRoot, "total time", all);
818}
819
820template <class T, typename DataType>
822 [[maybe_unused]] GeomMapView<T> &geomMap, [[maybe_unused]] int id,
823 [[maybe_unused]] DataType *data, [[maybe_unused]] CurveSharedPtr curve)
824{
825}
826
827template <>
829 [[maybe_unused]] GeomMapView<PointGeom> &geomMap, int id, NekDouble *data,
830 [[maybe_unused]] CurveSharedPtr curve)
831{
832 m_meshGraph->CreatePointGeom(m_meshGraph->GetSpaceDimension(), id, data[0],
833 data[1], data[2]);
834}
835
836template <>
838 [[maybe_unused]] GeomMapView<SegGeom> &geomMap, int id, int *data,
839 CurveSharedPtr curve)
840{
841 std::array<PointGeom *, 2> pts = {m_meshGraph->GetPointGeom(data[0]),
842 m_meshGraph->GetPointGeom(data[1])};
843
844 m_meshGraph->CreateSegGeom(id, m_meshGraph->GetSpaceDimension(), pts,
845 curve);
846}
847
848template <>
850 [[maybe_unused]] GeomMapView<TriGeom> &geomMap, int id, int *data,
851 CurveSharedPtr curve)
852{
853 std::array<SegGeom *, 3> segs = {m_meshGraph->GetSegGeom(data[0]),
854 m_meshGraph->GetSegGeom(data[1]),
855 m_meshGraph->GetSegGeom(data[2])};
856 m_meshGraph->CreateTriGeom(id, segs, curve);
857}
858
859template <>
861 [[maybe_unused]] GeomMapView<QuadGeom> &geomMap, int id, int *data,
862 CurveSharedPtr curve)
863{
864 std::array<SegGeom *, 4> segs = {
865 m_meshGraph->GetSegGeom(data[0]), m_meshGraph->GetSegGeom(data[1]),
866 m_meshGraph->GetSegGeom(data[2]), m_meshGraph->GetSegGeom(data[3])};
867 m_meshGraph->CreateQuadGeom(id, segs, curve);
868}
869
870template <>
872 [[maybe_unused]] GeomMapView<TetGeom> &geomMap, int id, int *data,
873 [[maybe_unused]] CurveSharedPtr curve)
874{
875 std::array<TriGeom *, 4> faces = {
876 m_meshGraph->GetTriGeom(data[0]), m_meshGraph->GetTriGeom(data[1]),
877 m_meshGraph->GetTriGeom(data[2]), m_meshGraph->GetTriGeom(data[3])};
878
879 auto geom = m_meshGraph->CreateTetGeom(id, faces);
880 m_meshGraph->PopulateFaceToElMap(geom, TetGeom::kNfaces);
881}
882
883template <>
885 [[maybe_unused]] GeomMapView<PyrGeom> &geomMap, int id, int *data,
886 [[maybe_unused]] CurveSharedPtr curve)
887{
888 std::array<Geometry2D *, 5> faces = {m_meshGraph->GetGeometry2D(data[0]),
889 m_meshGraph->GetGeometry2D(data[1]),
890 m_meshGraph->GetGeometry2D(data[2]),
891 m_meshGraph->GetGeometry2D(data[3]),
892 m_meshGraph->GetGeometry2D(data[4])};
893
894 auto geom = m_meshGraph->CreatePyrGeom(id, faces);
895 m_meshGraph->PopulateFaceToElMap(geom, PyrGeom::kNfaces);
896}
897
898template <>
900 [[maybe_unused]] GeomMapView<PrismGeom> &geomMap, int id, int *data,
901 [[maybe_unused]] CurveSharedPtr curve)
902{
903 std::array<Geometry2D *, 5> faces = {m_meshGraph->GetGeometry2D(data[0]),
904 m_meshGraph->GetGeometry2D(data[1]),
905 m_meshGraph->GetGeometry2D(data[2]),
906 m_meshGraph->GetGeometry2D(data[3]),
907 m_meshGraph->GetGeometry2D(data[4])};
908
909 auto geom = m_meshGraph->CreatePrismGeom(id, faces);
910 m_meshGraph->PopulateFaceToElMap(geom, PrismGeom::kNfaces);
911}
912
913template <>
915 [[maybe_unused]] GeomMapView<HexGeom> &geomMap, int id, int *data,
916 [[maybe_unused]] CurveSharedPtr curve)
917{
918 std::array<QuadGeom *, 6> faces = {
919 m_meshGraph->GetQuadGeom(data[0]), m_meshGraph->GetQuadGeom(data[1]),
920 m_meshGraph->GetQuadGeom(data[2]), m_meshGraph->GetQuadGeom(data[3]),
921 m_meshGraph->GetQuadGeom(data[4]), m_meshGraph->GetQuadGeom(data[5])};
922
923 auto geom = m_meshGraph->CreateHexGeom(id, faces);
924 m_meshGraph->PopulateFaceToElMap(geom, HexGeom::kNfaces);
925}
926
927template <class T, typename DataType>
929 const CurveMap &curveMap,
930 std::vector<int> &ids,
931 std::vector<DataType> &geomData)
932{
933 const int nGeomData = GetGeomDataDim(geomMap);
934 const int nRows = geomData.size() / nGeomData;
935 CurveSharedPtr empty;
936
937 // Construct geometry object.
938 if (curveMap.size() > 0)
939 {
940 for (int i = 0, cnt = 0; i < nRows; i++, cnt += nGeomData)
941 {
942 auto cIt = curveMap.find(ids[i]);
943 ConstructGeomObject(geomMap, ids[i], &geomData[cnt],
944 cIt == curveMap.end() ? empty : cIt->second);
945 }
946 }
947 else
948 {
949 for (int i = 0, cnt = 0; i < nRows; i++, cnt += nGeomData)
950 {
951 ConstructGeomObject(geomMap, ids[i], &geomData[cnt], empty);
952 }
953 }
954}
955
956template <class T, typename DataType>
958 std::string dataSet,
959 const std::unordered_set<int> &readIds,
960 std::vector<int> &ids,
961 std::vector<DataType> &geomData)
962{
963 if (!m_mesh->ContainsDataSet(dataSet))
964 {
965 return;
966 }
967
968 // Open mesh dataset
969 H5::DataSetSharedPtr data = m_mesh->OpenDataSet(dataSet);
970 H5::DataSpaceSharedPtr space = data->GetSpace();
971 std::vector<hsize_t> dims = space->GetDims();
972
973 // Open metadata dataset
974 H5::DataSetSharedPtr mdata = m_maps->OpenDataSet(dataSet);
975 H5::DataSpaceSharedPtr mspace = mdata->GetSpace();
976 std::vector<hsize_t> mdims = mspace->GetDims();
977
978 ASSERTL0(mdims[0] == dims[0], "map and data set lengths do not match");
979
980 const int nGeomData = GetGeomDataDim(geomMap);
981
982 // Read all IDs
983 std::vector<int> allIds;
984 mdata->Read(allIds, mspace);
985
986 // Selective reading; clear data space range so that we can select
987 // certain rows from the datasets.
988 space->ClearRange();
989
990 int i = 0;
991 std::vector<hsize_t> coords;
992 for (auto &id : allIds)
993 {
994 if (readIds.find(id) != readIds.end())
995 {
996 for (int j = 0; j < nGeomData; ++j)
997 {
998 coords.push_back(i);
999 coords.push_back(j);
1000 }
1001 ids.push_back(id);
1002 }
1003 ++i;
1004 }
1005
1006 space->SetSelection(coords.size() / 2, coords);
1007
1008 // Read selected data.
1009 data->Read(geomData, space, m_readPL);
1010}
1011
1012void MeshGraphIOHDF5::ReadCurveMap(CurveMap &curveMap, std::string dsName,
1013 const std::unordered_set<int> &readIds)
1014{
1015 auto &curveNodes = m_meshGraph->GetAllCurveNodes();
1016
1017 // If dataset does not exist, exit.
1018 if (!m_mesh->ContainsDataSet(dsName))
1019 {
1020 return;
1021 }
1022
1023 // Open up curve map data.
1024 H5::DataSetSharedPtr curveData = m_mesh->OpenDataSet(dsName);
1025 H5::DataSpaceSharedPtr curveSpace = curveData->GetSpace();
1026
1027 // Open up ID data set.
1028 H5::DataSetSharedPtr idData = m_maps->OpenDataSet(dsName);
1029 H5::DataSpaceSharedPtr idSpace = idData->GetSpace();
1030
1031 // Read all IDs and clear data space.
1032 std::vector<int> ids, newIds;
1033 idData->Read(ids, idSpace);
1034 curveSpace->ClearRange();
1035
1036 // Search IDs to figure out which curves to read.
1037 std::vector<hsize_t> curveSel;
1038
1039 int cnt = 0;
1040 for (auto &id : ids)
1041 {
1042 if (readIds.find(id) != readIds.end())
1043 {
1044 curveSel.push_back(cnt);
1045 curveSel.push_back(0);
1046 curveSel.push_back(cnt);
1047 curveSel.push_back(1);
1048 curveSel.push_back(cnt);
1049 curveSel.push_back(2);
1050 newIds.push_back(id);
1051 }
1052
1053 ++cnt;
1054 }
1055
1056 // Check to see whether any processor will read anything
1057 auto toRead = newIds.size();
1058 m_session->GetComm()->GetRowComm()->AllReduce(toRead,
1060
1061 if (toRead == 0)
1062 {
1063 return;
1064 }
1065
1066 // Now read curve map and read data.
1067 std::vector<int> curveInfo;
1068 curveSpace->SetSelection(curveSel.size() / 2, curveSel);
1069 curveData->Read(curveInfo, curveSpace, m_readPL);
1070
1071 curveSel.clear();
1072
1073 std::unordered_map<int, int> curvePtOffset;
1074
1075 // Construct curves. We'll populate nodes in a minute!
1076 for (int i = 0, cnt = 0, cnt2 = 0; i < curveInfo.size() / 3; ++i, cnt += 3)
1077 {
1079 newIds[i], (LibUtilities::PointsType)curveInfo[cnt + 1]);
1080
1081 curve->m_points.resize(curveInfo[cnt]);
1082
1083 const int ptOffset = curveInfo[cnt + 2];
1084
1085 for (int j = 0; j < curveInfo[cnt]; ++j)
1086 {
1087 // ptoffset gives us the row, multiply by 3 for number of
1088 // coordinates.
1089 curveSel.push_back(ptOffset + j);
1090 curveSel.push_back(0);
1091 curveSel.push_back(ptOffset + j);
1092 curveSel.push_back(1);
1093 curveSel.push_back(ptOffset + j);
1094 curveSel.push_back(2);
1095 }
1096
1097 // Store the offset so we know to come back later on to fill in
1098 // these points.
1099 curvePtOffset[newIds[i]] = 3 * cnt2;
1100 cnt2 += curveInfo[cnt];
1101
1102 curveMap[newIds[i]] = curve;
1103 }
1104
1105 curveInfo.clear();
1106
1107 // Open node data spacee.
1108 H5::DataSetSharedPtr nodeData = m_mesh->OpenDataSet("CURVE_NODES");
1109 H5::DataSpaceSharedPtr nodeSpace = nodeData->GetSpace();
1110
1111 nodeSpace->ClearRange();
1112 nodeSpace->SetSelection(curveSel.size() / 2, curveSel);
1113
1114 std::vector<NekDouble> nodeRawData;
1115 nodeData->Read(nodeRawData, nodeSpace, m_readPL);
1116
1117 // Go back and populate data from nodes.
1118 for (auto &cIt : curvePtOffset)
1119 {
1120 CurveSharedPtr curve = curveMap[cIt.first];
1121
1122 // Create nodes.
1123 int cnt = cIt.second;
1124 for (int i = 0; i < curve->m_points.size(); ++i, cnt += 3)
1125 {
1126 curveNodes.emplace_back(
1128 0, m_meshGraph->GetSpaceDimension(), nodeRawData[cnt],
1129 nodeRawData[cnt + 1], nodeRawData[cnt + 2]));
1130 curve->m_points[i] = curveNodes.back().get();
1131 }
1132 }
1133}
1134
1136{
1137 auto &domain = m_meshGraph->GetDomain();
1138
1139 if (m_inFormatVersion == 1)
1140 {
1141 std::map<int, CompositeSharedPtr> fullDomain;
1142 H5::DataSetSharedPtr dst = m_mesh->OpenDataSet("DOMAIN");
1143 H5::DataSpaceSharedPtr space = dst->GetSpace();
1144
1145 std::vector<std::string> data;
1146 dst->ReadVectorString(data, space, m_readPL);
1147 m_meshGraph->GetCompositeList(data[0], fullDomain);
1148 domain[0] = fullDomain;
1149
1150 return;
1151 }
1152
1153 std::vector<CompositeMap> fullDomain;
1154 H5::DataSetSharedPtr dst = m_mesh->OpenDataSet("DOMAIN");
1155 H5::DataSpaceSharedPtr space = dst->GetSpace();
1156
1157 std::vector<std::string> data;
1158 dst->ReadVectorString(data, space, m_readPL);
1159 for (auto &dIt : data)
1160 {
1161 fullDomain.push_back(CompositeMap());
1162 m_meshGraph->GetCompositeList(dIt, fullDomain.back());
1163 }
1164
1165 H5::DataSetSharedPtr mdata = m_maps->OpenDataSet("DOMAIN");
1166 H5::DataSpaceSharedPtr mspace = mdata->GetSpace();
1167
1168 std::vector<int> ids;
1169 mdata->Read(ids, mspace);
1170
1171 for (int i = 0; i < ids.size(); ++i)
1172 {
1173 domain[ids[i]] = fullDomain[i];
1174 }
1175}
1176
1178{
1179 if (!rng || rng->m_compElmts == false)
1180 {
1181 return; // composite range not being used.
1182 }
1183
1184 std::string nm = "COMPOSITE";
1185
1186 H5::DataSetSharedPtr data = m_mesh->OpenDataSet(nm);
1187 H5::DataSpaceSharedPtr space = data->GetSpace();
1188 std::vector<hsize_t> dims = space->GetDims();
1189
1190 std::vector<std::string> comps;
1191 data->ReadVectorString(comps, space);
1192
1193 H5::DataSetSharedPtr mdata = m_maps->OpenDataSet(nm);
1194 H5::DataSpaceSharedPtr mspace = mdata->GetSpace();
1195 std::vector<hsize_t> mdims = mspace->GetDims();
1196
1197 std::vector<int> ids;
1198 mdata->Read(ids, mspace);
1199
1200 for (int i = 0; i < dims[0]; i++)
1201 {
1202
1203 if (rng->m_comps.count(ids[i]))
1204 {
1205
1206 std::string compStr = comps[i];
1207
1208 char type;
1209 std::istringstream strm(compStr);
1210
1211 strm >> type;
1212
1213 CompositeSharedPtr comp =
1215
1216 std::string::size_type indxBeg = compStr.find_first_of('[') + 1;
1217 std::string::size_type indxEnd = compStr.find_last_of(']') - 1;
1218
1219 std::string indxStr =
1220 compStr.substr(indxBeg, indxEnd - indxBeg + 1);
1221 std::vector<unsigned int> seqVector;
1222
1223 ParseUtils::GenerateSeqVector(indxStr, seqVector);
1224
1225 for (auto it : seqVector) // add ids to Traceid
1226 {
1227 rng->m_traceIDs.insert(it);
1228 }
1229 }
1230 }
1231}
1232
1234{
1235 auto &vertSet = m_meshGraph->GetGeomMap<PointGeom>();
1236 auto &segGeoms = m_meshGraph->GetGeomMap<SegGeom>();
1237 auto &triGeoms = m_meshGraph->GetGeomMap<TriGeom>();
1238 auto &quadGeoms = m_meshGraph->GetGeomMap<QuadGeom>();
1239 auto &tetGeoms = m_meshGraph->GetGeomMap<TetGeom>();
1240 auto &pyrGeoms = m_meshGraph->GetGeomMap<PyrGeom>();
1241 auto &prismGeoms = m_meshGraph->GetGeomMap<PrismGeom>();
1242 auto &hexGeoms = m_meshGraph->GetGeomMap<HexGeom>();
1243 CompositeMap &meshComposites = m_meshGraph->GetComposites();
1244
1245 std::string nm = "COMPOSITE";
1246
1247 H5::DataSetSharedPtr data = m_mesh->OpenDataSet(nm);
1248 H5::DataSpaceSharedPtr space = data->GetSpace();
1249 std::vector<hsize_t> dims = space->GetDims();
1250
1251 std::vector<std::string> comps;
1252 data->ReadVectorString(comps, space);
1253
1254 H5::DataSetSharedPtr mdata = m_maps->OpenDataSet(nm);
1255 H5::DataSpaceSharedPtr mspace = mdata->GetSpace();
1256 std::vector<hsize_t> mdims = mspace->GetDims();
1257
1258 std::vector<int> ids;
1259 mdata->Read(ids, mspace);
1260
1261 auto &graph_compOrder = m_meshGraph->GetCompositeOrdering();
1262 for (int i = 0; i < dims[0]; i++)
1263 {
1264 std::string compStr = comps[i];
1265
1266 char type;
1267 std::istringstream strm(compStr);
1268
1269 strm >> type;
1270
1272
1273 std::string::size_type indxBeg = compStr.find_first_of('[') + 1;
1274 std::string::size_type indxEnd = compStr.find_last_of(']') - 1;
1275
1276 std::string indxStr = compStr.substr(indxBeg, indxEnd - indxBeg + 1);
1277 std::vector<unsigned int> seqVector;
1278
1279 ParseUtils::GenerateSeqVector(indxStr, seqVector);
1280 m_compOrder[ids[i]] = seqVector;
1281 graph_compOrder[ids[i]] = seqVector;
1282
1283 switch (type)
1284 {
1285 case 'V':
1286 for (auto &i : seqVector)
1287 {
1288 auto it = vertSet.find(i);
1289 if (it != vertSet.end())
1290 {
1291 comp->m_geomVec.push_back((*it).second);
1292 }
1293 }
1294 break;
1295 case 'S':
1296 case 'E':
1297 for (auto &i : seqVector)
1298 {
1299 auto it = segGeoms.find(i);
1300 if (it != segGeoms.end())
1301 {
1302 comp->m_geomVec.push_back((*it).second);
1303 }
1304 }
1305 break;
1306 case 'Q':
1307 for (auto &i : seqVector)
1308 {
1309 auto it = quadGeoms.find(i);
1310 if (it != quadGeoms.end())
1311 {
1312 if (m_meshGraph->CheckRange(*(*it).second))
1313 {
1314 comp->m_geomVec.push_back((*it).second);
1315 }
1316 }
1317 }
1318 break;
1319 case 'T':
1320 for (auto &i : seqVector)
1321 {
1322 auto it = triGeoms.find(i);
1323 if (it != triGeoms.end())
1324 {
1325 if (m_meshGraph->CheckRange(*(*it).second))
1326 {
1327 comp->m_geomVec.push_back((*it).second);
1328 }
1329 }
1330 }
1331 break;
1332 case 'F':
1333 for (auto &i : seqVector)
1334 {
1335 auto it1 = quadGeoms.find(i);
1336 if (it1 != quadGeoms.end())
1337 {
1338 if (m_meshGraph->CheckRange(*(*it1).second))
1339 {
1340 comp->m_geomVec.push_back((*it1).second);
1341 }
1342 }
1343 auto it2 = triGeoms.find(i);
1344 if (it2 != triGeoms.end())
1345 {
1346 if (m_meshGraph->CheckRange(*(*it2).second))
1347 {
1348 comp->m_geomVec.push_back((*it2).second);
1349 }
1350 }
1351 }
1352 break;
1353 case 'A':
1354 for (auto &i : seqVector)
1355 {
1356 auto it = tetGeoms.find(i);
1357 if (it != tetGeoms.end())
1358 {
1359 if (m_meshGraph->CheckRange(*(*it).second))
1360 {
1361 comp->m_geomVec.push_back((*it).second);
1362 }
1363 }
1364 }
1365 break;
1366 case 'P':
1367 for (auto &i : seqVector)
1368 {
1369 auto it = pyrGeoms.find(i);
1370 if (it != pyrGeoms.end())
1371 {
1372 if (m_meshGraph->CheckRange(*(*it).second))
1373 {
1374 comp->m_geomVec.push_back((*it).second);
1375 }
1376 }
1377 }
1378 break;
1379 case 'R':
1380 for (auto &i : seqVector)
1381 {
1382 auto it = prismGeoms.find(i);
1383 if (it != prismGeoms.end())
1384 {
1385 if (m_meshGraph->CheckRange(*(*it).second))
1386 {
1387 comp->m_geomVec.push_back((*it).second);
1388 }
1389 }
1390 }
1391 break;
1392 case 'H':
1393 for (auto &i : seqVector)
1394 {
1395 auto it = hexGeoms.find(i);
1396 if (it != hexGeoms.end())
1397 {
1398 if (m_meshGraph->CheckRange(*(*it).second))
1399 {
1400 comp->m_geomVec.push_back((*it).second);
1401 }
1402 }
1403 }
1404 break;
1405 }
1406
1407 if (comp->m_geomVec.size() > 0)
1408 {
1409 meshComposites[ids[i]] = comp;
1410 }
1411 }
1412}
1413
1415 std::unordered_map<int, int> &id2row)
1416{
1418
1419 std::string nm = "COMPOSITE";
1420
1421 H5::DataSetSharedPtr data = m_mesh->OpenDataSet(nm);
1422 H5::DataSpaceSharedPtr space = data->GetSpace();
1423 std::vector<hsize_t> dims = space->GetDims();
1424
1425 std::vector<std::string> comps;
1426 data->ReadVectorString(comps, space);
1427
1428 H5::DataSetSharedPtr mdata = m_maps->OpenDataSet(nm);
1429 H5::DataSpaceSharedPtr mspace = mdata->GetSpace();
1430 std::vector<hsize_t> mdims = mspace->GetDims();
1431
1432 std::vector<int> ids;
1433 mdata->Read(ids, mspace);
1434
1435 for (int i = 0; i < dims[0]; i++)
1436 {
1437 std::string compStr = comps[i];
1438
1439 char type;
1440 std::istringstream strm(compStr);
1441
1442 strm >> type;
1443
1444 std::string::size_type indxBeg = compStr.find_first_of('[') + 1;
1445 std::string::size_type indxEnd = compStr.find_last_of(']') - 1;
1446
1447 std::string indxStr = compStr.substr(indxBeg, indxEnd - indxBeg + 1);
1448 std::vector<unsigned int> seqVector;
1449 ParseUtils::GenerateSeqVector(indxStr, seqVector);
1450
1452
1453 switch (type)
1454 {
1455 case 'V':
1456 shapeType = LibUtilities::ePoint;
1457 break;
1458 case 'S':
1459 case 'E':
1460 shapeType = LibUtilities::eSegment;
1461 break;
1462 case 'Q':
1463 case 'F':
1464 // Note that for HDF5, the composite descriptor is only used
1465 // for partitioning purposes so 'F' tag is not really going
1466 // to be critical in this context.
1467 shapeType = LibUtilities::eQuadrilateral;
1468 break;
1469 case 'T':
1470 shapeType = LibUtilities::eTriangle;
1471 break;
1472 case 'A':
1473 shapeType = LibUtilities::eTetrahedron;
1474 break;
1475 case 'P':
1476 shapeType = LibUtilities::ePyramid;
1477 break;
1478 case 'R':
1479 shapeType = LibUtilities::ePrism;
1480 break;
1481 case 'H':
1482 shapeType = LibUtilities::eHexahedron;
1483 break;
1484 }
1485
1486 ASSERTL0(shapeType != eNoShapeType, "Invalid shape.");
1487
1488 std::vector<int> filteredVector;
1489 for (auto &compElmt : seqVector)
1490 {
1491 if (id2row.find(compElmt) == id2row.end())
1492 {
1493 continue;
1494 }
1495
1496 filteredVector.push_back(compElmt);
1497 }
1498
1499 if (filteredVector.size() == 0)
1500 {
1501 continue;
1502 }
1503
1504 ret[ids[i]] = std::make_pair(shapeType, filteredVector);
1505 }
1506
1507 return ret;
1508}
1509
1510template <class T, typename std::enable_if<T::kDim == 0, int>::type = 0>
1511inline NekDouble GetGeomData(T *geom, int i)
1512{
1513 return (*geom)(i);
1514}
1515
1516template <class T, typename std::enable_if<T::kDim == 1, int>::type = 0>
1517inline int GetGeomData(T *geom, int i)
1518{
1519 return geom->GetVid(i);
1520}
1521
1522template <class T, typename std::enable_if<T::kDim == 2, int>::type = 0>
1523inline int GetGeomData(T *geom, int i)
1524{
1525 return geom->GetEid(i);
1526}
1527
1528template <class T, typename std::enable_if<T::kDim == 3, int>::type = 0>
1529inline int GetGeomData(T *geom, int i)
1530{
1531 return geom->GetFid(i);
1532}
1533
1534template <class T>
1536 std::string datasetName)
1537{
1538 typedef typename std::conditional<std::is_same_v<T, PointGeom>, NekDouble,
1539 int>::type DataType;
1540
1541 const int nGeomData = GetGeomDataDim(geomMap);
1542 const size_t nGeom = geomMap.size();
1543
1544 if (nGeom == 0)
1545 {
1546 return;
1547 }
1548
1549 // Construct a map storing IDs
1550 std::vector<int> idMap(nGeom);
1551 std::vector<DataType> data(nGeom * nGeomData);
1552
1553 int cnt1 = 0, cnt2 = 0;
1554 for (auto [id, geom] : geomMap)
1555 {
1556 idMap[cnt1++] = id;
1557
1558 for (int j = 0; j < nGeomData; ++j)
1559 {
1560 data[cnt2 + j] = GetGeomData(geom, j);
1561 }
1562
1563 cnt2 += nGeomData;
1564 }
1565
1566 std::vector<hsize_t> dims = {static_cast<hsize_t>(nGeom),
1567 static_cast<hsize_t>(nGeomData)};
1570 std::shared_ptr<H5::DataSpace>(new H5::DataSpace(dims));
1571 H5::DataSetSharedPtr dst = m_mesh->CreateDataSet(datasetName, tp, ds);
1572 dst->Write(data, ds);
1573
1574 tp = H5::DataType::OfObject(idMap[0]);
1575 dims = {nGeom};
1576 ds = std::shared_ptr<H5::DataSpace>(new H5::DataSpace(dims));
1577 dst = m_maps->CreateDataSet(datasetName, tp, ds);
1578 dst->Write(idMap, ds);
1579}
1580
1581void MeshGraphIOHDF5::WriteCurveMap(CurveMap &curves, std::string dsName,
1582 MeshCurvedPts &curvedPts, int &ptOffset,
1583 int &newIdx)
1584{
1585 std::vector<int> data, map;
1586
1587 // Compile curve data.
1588 for (auto &c : curves)
1589 {
1590 map.push_back(c.first);
1591 data.push_back(c.second->m_points.size());
1592 data.push_back(c.second->m_ptype);
1593 data.push_back(ptOffset);
1594
1595 ptOffset += c.second->m_points.size();
1596
1597 for (auto &pt : c.second->m_points)
1598 {
1599 MeshVertex v;
1600 v.id = newIdx;
1601 pt->GetCoords(v.x, v.y, v.z);
1602 curvedPts.pts.push_back(v);
1603 curvedPts.index.push_back(newIdx++);
1604 }
1605 }
1606
1607 // Write data.
1608 std::vector<hsize_t> dims = {data.size() / 3, 3};
1611 std::shared_ptr<H5::DataSpace>(new H5::DataSpace(dims));
1612 H5::DataSetSharedPtr dst = m_mesh->CreateDataSet(dsName, tp, ds);
1613 dst->Write(data, ds);
1614
1615 tp = H5::DataType::OfObject(map[0]);
1616 dims = {map.size()};
1617 ds = std::shared_ptr<H5::DataSpace>(new H5::DataSpace(dims));
1618 dst = m_maps->CreateDataSet(dsName, tp, ds);
1619 dst->Write(map, ds);
1620}
1621
1623{
1624 std::vector<double> vertData(curvedPts.pts.size() * 3);
1625
1626 int cnt = 0;
1627 for (auto &pt : curvedPts.pts)
1628 {
1629 vertData[cnt++] = pt.x;
1630 vertData[cnt++] = pt.y;
1631 vertData[cnt++] = pt.z;
1632 }
1633
1634 std::vector<hsize_t> dims = {curvedPts.pts.size(), 3};
1637 std::shared_ptr<H5::DataSpace>(new H5::DataSpace(dims));
1638 H5::DataSetSharedPtr dst = m_mesh->CreateDataSet("CURVE_NODES", tp, ds);
1639 dst->Write(vertData, ds);
1640}
1641
1643{
1644 std::vector<std::string> comps;
1645
1646 // dont need location map only a id map
1647 // will filter the composites per parition on read, its easier
1648 // composites do not need to be written in paralell.
1649 std::vector<int> c_map;
1650
1651 for (auto &cIt : composites)
1652 {
1653 if (cIt.second->m_geomVec.size() == 0)
1654 {
1655 continue;
1656 }
1657
1658 comps.push_back(GetCompositeString(cIt.second));
1659 c_map.push_back(cIt.first);
1660 }
1661
1664 H5::DataSetSharedPtr dst = m_mesh->CreateDataSet("COMPOSITE", tp, ds);
1665 dst->WriteVectorString(comps, ds, tp);
1666
1667 tp = H5::DataType::OfObject(c_map[0]);
1668 ds = H5::DataSpace::OneD(c_map.size());
1669 dst = m_maps->CreateDataSet("COMPOSITE", tp, ds);
1670 dst->Write(c_map, ds);
1671}
1672
1673void MeshGraphIOHDF5::WriteDomain(std::map<int, CompositeMap> &domain)
1674{
1675 // dont need location map only a id map
1676 // will filter the composites per parition on read, its easier
1677 // composites do not need to be written in paralell.
1678 std::vector<int> d_map;
1679 std::vector<std::vector<unsigned int>> idxList;
1680
1681 int cnt = 0;
1682 for (auto &dIt : domain)
1683 {
1684 idxList.push_back(std::vector<unsigned int>());
1685 for (auto cIt = dIt.second.begin(); cIt != dIt.second.end(); ++cIt)
1686 {
1687 idxList[cnt].push_back(cIt->first);
1688 }
1689
1690 ++cnt;
1691 d_map.push_back(dIt.first);
1692 }
1693
1694 std::stringstream domString;
1695 std::vector<std::string> doms;
1696 for (auto &cIt : idxList)
1697 {
1698 doms.push_back(ParseUtils::GenerateSeqString(cIt));
1699 }
1700
1703 H5::DataSetSharedPtr dst = m_mesh->CreateDataSet("DOMAIN", tp, ds);
1704 dst->WriteVectorString(doms, ds, tp);
1705
1706 tp = H5::DataType::OfObject(d_map[0]);
1707 ds = H5::DataSpace::OneD(d_map.size());
1708 dst = m_maps->CreateDataSet("DOMAIN", tp, ds);
1709 dst->Write(d_map, ds);
1710}
1711
1713 const std::string &outfilename, bool defaultExp,
1714 [[maybe_unused]] const LibUtilities::FieldMetaDataMap &metadata)
1715{
1716 auto &vertSet = m_meshGraph->GetGeomMap<PointGeom>();
1717 auto &segGeoms = m_meshGraph->GetGeomMap<SegGeom>();
1718 auto &triGeoms = m_meshGraph->GetGeomMap<TriGeom>();
1719 auto &quadGeoms = m_meshGraph->GetGeomMap<QuadGeom>();
1720 auto &hexGeoms = m_meshGraph->GetGeomMap<HexGeom>();
1721 auto &pyrGeoms = m_meshGraph->GetGeomMap<PyrGeom>();
1722 auto &prismGeoms = m_meshGraph->GetGeomMap<PrismGeom>();
1723 auto &tetGeoms = m_meshGraph->GetGeomMap<TetGeom>();
1724 CurveMap &curvedEdges = m_meshGraph->GetCurvedEdges();
1725 CurveMap &curvedFaces = m_meshGraph->GetCurvedFaces();
1726 CompositeMap &meshComposites = m_meshGraph->GetComposites();
1727 auto domain = m_meshGraph->GetDomain();
1728
1729 std::vector<std::string> tmp;
1730 boost::split(tmp, outfilename, boost::is_any_of("."));
1731 std::string filenameXml = tmp[0] + ".xml";
1732 std::string filenameHdf5 = tmp[0] + ".nekg";
1733
1734 //////////////////
1735 // XML part
1736 //////////////////
1737
1738 // Check to see if a xml of the same name exists
1739 // if might have boundary conditions etc, we will just alter the
1740 // geometry tag if needed
1741 TiXmlDocument *doc = new TiXmlDocument;
1742 TiXmlElement *root;
1743 TiXmlElement *geomTag;
1744
1745 if (fs::exists(filenameXml.c_str()))
1746 {
1747 std::ifstream file(filenameXml.c_str());
1748 file >> (*doc);
1749 TiXmlHandle docHandle(doc);
1750 root = docHandle.FirstChildElement("NEKTAR").Element();
1751 ASSERTL0(root, "Unable to find NEKTAR tag in file.");
1752 geomTag = root->FirstChildElement("GEOMETRY");
1753 defaultExp = false;
1754 }
1755 else
1756 {
1757 TiXmlDeclaration *decl = new TiXmlDeclaration("1.0", "utf-8", "");
1758 doc->LinkEndChild(decl);
1759 root = new TiXmlElement("NEKTAR");
1760 doc->LinkEndChild(root);
1761
1762 geomTag = new TiXmlElement("GEOMETRY");
1763 root->LinkEndChild(geomTag);
1764 }
1765
1766 // Update attributes with dimensions.
1767 geomTag->SetAttribute("DIM", m_meshGraph->GetMeshDimension());
1768 geomTag->SetAttribute("SPACE", m_meshGraph->GetSpaceDimension());
1769 geomTag->SetAttribute("HDF5FILE", filenameHdf5);
1770
1771 geomTag->Clear();
1772
1773 if (defaultExp)
1774 {
1775 TiXmlElement *expTag = new TiXmlElement("EXPANSIONS");
1776
1777 for (auto it = meshComposites.begin(); it != meshComposites.end(); it++)
1778 {
1779 if (it->second->m_geomVec[0]->GetShapeDim() ==
1780 m_meshGraph->GetMeshDimension())
1781 {
1782 TiXmlElement *exp = new TiXmlElement("E");
1783 exp->SetAttribute("COMPOSITE",
1784 "C[" + std::to_string(it->first) + "]");
1785 exp->SetAttribute("NUMMODES", 4);
1786 exp->SetAttribute("TYPE", "MODIFIED");
1787 exp->SetAttribute("FIELDS", "u");
1788
1789 expTag->LinkEndChild(exp);
1790 }
1791 }
1792 root->LinkEndChild(expTag);
1793 }
1794
1795 auto movement = m_meshGraph->GetMovement();
1796 if (movement)
1797 {
1798 movement->WriteMovement(root);
1799 }
1800
1801 doc->SaveFile(filenameXml);
1802
1803 //////////////////
1804 // HDF5 part
1805 //////////////////
1806
1807 // This is serial IO so we will just override any existing file.
1808 m_file = H5::File::Create(filenameHdf5, H5F_ACC_TRUNC);
1809 auto hdfRoot = m_file->CreateGroup("NEKTAR");
1810 auto hdfRoot2 = hdfRoot->CreateGroup("GEOMETRY");
1811
1812 // Write format version.
1813 hdfRoot2->SetAttribute("FORMAT_VERSION", FORMAT_VERSION);
1814
1815 // Create main groups.
1816 m_mesh = hdfRoot2->CreateGroup("MESH");
1817 m_maps = hdfRoot2->CreateGroup("MAPS");
1818
1819 WriteGeometryMap(vertSet, "VERT");
1820 WriteGeometryMap(segGeoms, "SEG");
1821 if (m_meshGraph->GetMeshDimension() > 1)
1822 {
1823 WriteGeometryMap(triGeoms, "TRI");
1824 WriteGeometryMap(quadGeoms, "QUAD");
1825 }
1826 if (m_meshGraph->GetMeshDimension() > 2)
1827 {
1828 WriteGeometryMap(tetGeoms, "TET");
1829 WriteGeometryMap(pyrGeoms, "PYR");
1830 WriteGeometryMap(prismGeoms, "PRISM");
1831 WriteGeometryMap(hexGeoms, "HEX");
1832 }
1833
1834 // Write curves
1835 int ptOffset = 0, newIdx = 0;
1836 MeshCurvedPts curvePts;
1837 WriteCurveMap(curvedEdges, "CURVE_EDGE", curvePts, ptOffset, newIdx);
1838 WriteCurveMap(curvedFaces, "CURVE_FACE", curvePts, ptOffset, newIdx);
1839 WriteCurvePoints(curvePts);
1840
1841 // Write composites and domain.
1842 WriteComposites(meshComposites);
1843 WriteDomain(domain);
1844}
1845
1846} // namespace Nektar::SpatialDomains
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define TIME_RESULT(verb, msg, timer)
HDF5 DataSpace wrapper.
Definition H5.h:309
static DataSpaceSharedPtr OneD(hsize_t size)
Definition H5.cpp:411
static DataTypeSharedPtr OfObject(const T &obj)
Definition H5.h:397
static DataTypeSharedPtr String(size_t len=0)
Definition H5.cpp:517
static FileSharedPtr Open(const std::string &filename, unsigned mode, PListSharedPtr accessPL=PList::Default())
Definition H5.cpp:636
static FileSharedPtr Create(const std::string &filename, unsigned mode, PListSharedPtr createPL=PList::Default(), PListSharedPtr accessPL=PList::Default())
Definition H5.cpp:627
static PListSharedPtr DatasetXfer()
Properties for raw data transfer.
Definition H5.cpp:109
static PListSharedPtr FileAccess()
Properties for file access.
Definition H5.cpp:91
static PListSharedPtr Default()
Default options.
Definition H5.cpp:73
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
bool ModuleExists(tKey idKey)
Checks if a particular module is available.
static void GetXMLElementTimeLevel(TiXmlElement *&element, const size_t timeLevel, const bool enableCheck=true)
Get XML elment time level (Parallel-in-Time)
static std::string RegisterCmdLineFlag(const std::string &pName, const std::string &pShortName, const std::string &pDescription)
Registers a command-line flag with the session reader.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Generic object pool allocator/deallocator.
static std::string GenerateSeqString(const std::vector< T > &v)
Generate a compressed comma-separated string representation of a vector of unsigned integers.
Definition ParseUtils.h:72
static bool GenerateSeqVector(const std::string &str, std::vector< unsigned int > &out)
Takes a comma-separated compressed string and converts it to entries in a vector.
static const int kNfaces
Definition HexGeom.h:52
void FillGeomMap(GeomMapView< T > &geomMap, const CurveMap &curveMap, std::vector< int > &ids, std::vector< DataType > &geomData)
void v_PartitionMesh(LibUtilities::SessionReaderSharedPtr session) final
Partition the mesh.
void SetupCompositeRange(LibUtilities::DomainRangeShPtr &rng)
LibUtilities::H5::GroupSharedPtr m_mesh
LibUtilities::H5::GroupSharedPtr m_maps
void ConstructGeomObject(GeomMapView< T > &geomMap, int id, DataType *data, CurveSharedPtr curve)
void v_ReadGeometry(bool fillGraph) final
void WriteDomain(std::map< int, CompositeMap > &domain)
static const unsigned int FORMAT_VERSION
Version of the Nektar++ HDF5 geometry format, which is embedded into the main NEKTAR/GEOMETRY group a...
void WriteCurveMap(CurveMap &curves, std::string dsName, MeshCurvedPts &curvedPts, int &ptOffset, int &newIdx)
void v_WriteGeometry(const std::string &outfilename, bool defaultExp=false, const LibUtilities::FieldMetaDataMap &metadata=LibUtilities::NullFieldMetaDataMap) final
void WriteCurvePoints(MeshCurvedPts &curvedPts)
LibUtilities::H5::FileSharedPtr m_file
void ReadCurveMap(CurveMap &curveMap, std::string dsName, const std::unordered_set< int > &readIds)
void ReadGeometryData(GeomMapView< T > &geomMap, std::string dataSet, const std::unordered_set< int > &readIds, std::vector< int > &ids, std::vector< DataType > &geomData)
void WriteGeometryMap(GeomMapView< T > &geomMap, std::string datasetName)
LibUtilities::H5::PListSharedPtr m_readPL
static MeshGraphIOSharedPtr create()
LibUtilities::SessionReaderSharedPtr m_session
Definition MeshGraphIO.h:88
std::string GetCompositeString(CompositeSharedPtr comp)
Returns a string representation of a composite.
CompositeDescriptor CreateCompositeDescriptor()
static const int kNfaces
Definition PyrGeom.h:52
static const int kNfaces
Definition TetGeom.h:52
std::shared_ptr< DataSpace > DataSpaceSharedPtr
Definition H5.h:78
std::shared_ptr< PList > PListSharedPtr
Definition H5.h:92
std::shared_ptr< File > FileSharedPtr
Definition H5.h:88
std::shared_ptr< DataType > DataTypeSharedPtr
Definition H5.h:74
std::shared_ptr< Group > GroupSharedPtr
Definition FieldIOHdf5.h:48
std::shared_ptr< DataSet > DataSetSharedPtr
Definition H5.h:90
std::map< std::string, std::string > FieldMetaDataMap
Definition FieldIO.h:50
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< DomainRange > DomainRangeShPtr
Definition DomainRange.h:69
static DomainRangeShPtr NullDomainRangeShPtr
Definition DomainRange.h:70
CommFactory & GetCommFactory()
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition Comm.h:55
NekDouble GetGeomData(T *geom, int i)
std::map< int, std::pair< LibUtilities::ShapeType, std::vector< int > > > CompositeDescriptor
Definition MeshGraph.h:124
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< Composite > CompositeSharedPtr
Definition MeshGraph.h:178
std::shared_ptr< Curve > CurveSharedPtr
Definition Curve.hpp:60
std::unordered_map< int, CurveSharedPtr > CurveMap
Definition Curve.hpp:61
MeshGraphIOFactory & GetMeshGraphIOFactory()
int GetGeomDataDim(GeomMapView< T > &geomMap)
std::shared_ptr< MeshPartition > MeshPartitionSharedPtr
void UniqueValues(std::unordered_set< int > &unique)
std::map< int, CompositeSharedPtr > CompositeMap
Definition MeshGraph.h:179
std::vector< int64_t > index
Mapping to access the pts value. Given a 'ptoffset' value the npoints subsquent values provide the in...
std::vector< MeshVertex > pts
mapping to access pts value.
std::vector< unsigned int > list