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