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