Nektar++
MeshGraphIOXml.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: MeshGraphIOXml.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:
32//
33//
34////////////////////////////////////////////////////////////////////////////////
35
36#include <iomanip>
37
41
46
47#include <boost/format.hpp>
48
49#include <tinyxml.h>
50
51using namespace std;
52
54{
55
56std::string MeshGraphIOXml::className =
58 "Xml", MeshGraphIOXml::create, "IO with Xml geometry");
59
62{
63 // Get row of comm, or the whole comm if not split
64 LibUtilities::CommSharedPtr comm = session->GetComm();
65 LibUtilities::CommSharedPtr commMesh = comm->GetRowComm();
66 const bool isRoot = comm->TreatAsRankZero();
67
68 m_session = session;
69
70 int meshDimension;
71
72 // Load file for root process only (since this is always needed)
73 // and determine if the provided geometry has already been
74 // partitioned. This will be the case if the user provides the
75 // directory of mesh partitions as an input. Partitioned geometries
76 // have the attribute
77 // PARTITION=X
78 // where X is the number of the partition (and should match the
79 // process rank). The result is shared with all other processes.
80 int isPartitioned = 0;
81 if (isRoot)
82 {
83 if (m_session->DefinesElement("Nektar/Geometry"))
84 {
85 if (m_session->GetElement("Nektar/Geometry")
86 ->Attribute("PARTITION"))
87 {
88 std::cout << "Using pre-partitioned mesh." << std::endl;
89 isPartitioned = 1;
90 }
91 m_session->GetElement("NEKTAR/GEOMETRY")
92 ->QueryIntAttribute("DIM", &meshDimension);
93 }
94 }
95 comm->Bcast(isPartitioned, 0);
96 comm->Bcast(meshDimension, 0);
97
98 // If the mesh is already partitioned, we are done. Remaining
99 // processes must load their partitions.
100 if (isPartitioned)
101 {
102 if (!isRoot)
103 {
104 m_session->InitSession();
105 }
106 }
107 else
108 {
109 // Default partitioner to use is Scotch, if it is installed.
110 bool haveScotch = GetMeshPartitionFactory().ModuleExists("Scotch");
111 bool haveMetis = GetMeshPartitionFactory().ModuleExists("Metis");
112
113 string partitionerName = "Scotch";
114 if (!haveScotch && haveMetis)
115 {
116 partitionerName = "Metis";
117 }
118
119 // Override default with command-line flags if they are set.
120 if (session->DefinesCmdLineArgument("use-metis"))
121 {
122 partitionerName = "Metis";
123 }
124 if (session->DefinesCmdLineArgument("use-scotch"))
125 {
126 partitionerName = "Scotch";
127 }
128
129 // Mesh has not been partitioned so do partitioning if required. Note
130 // in the serial case nothing is done as we have already loaded the
131 // mesh.
132 if (session->DefinesCmdLineArgument("part-only") ||
133 session->DefinesCmdLineArgument("part-only-overlapping"))
134 {
135 // Perform partitioning of the mesh only. For this we insist the
136 // code is run in serial (parallel execution is pointless).
137 ASSERTL0(comm->GetSize() == 1,
138 "The 'part-only' option should be used in serial.");
139
140 // Read 'lite' geometry information
142
143 // Number of partitions is specified by the parameter.
144 int nParts;
145 auto comp = CreateCompositeDescriptor();
146
147 MeshPartitionSharedPtr partitioner =
149 partitionerName, session, session->GetComm(), meshDimension,
150 CreateMeshEntities(), comp);
151
152 if (session->DefinesCmdLineArgument("part-only"))
153 {
154 nParts = session->GetCmdLineArgument<int>("part-only");
155 partitioner->PartitionMesh(nParts, true);
156 }
157 else
158 {
159 nParts =
160 session->GetCmdLineArgument<int>("part-only-overlapping");
161 partitioner->PartitionMesh(nParts, true, true);
162 }
163
164 vector<set<unsigned int>> elmtIDs;
165 vector<unsigned int> parts(nParts);
166 for (int i = 0; i < nParts; ++i)
167 {
168 vector<unsigned int> elIDs;
169 set<unsigned int> tmp;
170 partitioner->GetElementIDs(i, elIDs);
171 tmp.insert(elIDs.begin(), elIDs.end());
172 elmtIDs.push_back(tmp);
173 parts[i] = i;
174 }
175
176 this->WriteXMLGeometry(m_session->GetSessionName(), elmtIDs, parts);
177
178 if (isRoot && session->DefinesCmdLineArgument("part-info"))
179 {
180 partitioner->PrintPartInfo(std::cout);
181 }
182
183 session->Finalise();
184 exit(0);
185 }
186
187 if (commMesh->GetSize() > 1)
188 {
189 ASSERTL0(haveScotch || haveMetis,
190 "Valid partitioner not found! Either Scotch or METIS "
191 "should be used.");
192
193 int nParts = commMesh->GetSize();
194
195 if (session->GetSharedFilesystem())
196 {
197 vector<unsigned int> keys, vals;
198 int i;
199
200 if (isRoot)
201 {
202 // Read 'lite' geometry information
204
205 // Store composite ordering and boundary information.
207 m_meshGraph->SetCompositeOrdering(m_compOrder);
208 auto comp = CreateCompositeDescriptor();
209
210 // Create mesh partitioner.
211 MeshPartitionSharedPtr partitioner =
213 partitionerName, session, session->GetComm(),
214 meshDimension, CreateMeshEntities(), comp);
215
216 partitioner->PartitionMesh(nParts, true);
217
218 vector<set<unsigned int>> elmtIDs;
219 vector<unsigned int> parts(nParts);
220 for (i = 0; i < nParts; ++i)
221 {
222 vector<unsigned int> elIDs;
223 set<unsigned int> tmp;
224 partitioner->GetElementIDs(i, elIDs);
225 tmp.insert(elIDs.begin(), elIDs.end());
226 elmtIDs.push_back(tmp);
227 parts[i] = i;
228 }
229
230 // Call WriteGeometry to write out partition files. This
231 // will populate m_bndRegOrder.
232 this->WriteXMLGeometry(m_session->GetSessionName(), elmtIDs,
233 parts);
234
235 // Communicate orderings to the other processors.
236
237 // First send sizes of the orderings and boundary
238 // regions to allocate storage on the remote end.
239 keys.resize(2);
240 keys[0] = m_compOrder.size();
241 keys[1] = m_bndRegOrder.size();
242 comm->Bcast(keys, 0);
243
244 // Construct the keys and sizes of values for composite
245 // ordering
246 keys.resize(m_compOrder.size());
247 vals.resize(m_compOrder.size());
248
249 i = 0;
250 for (auto &cIt : m_compOrder)
251 {
252 keys[i] = cIt.first;
253 vals[i++] = cIt.second.size();
254 }
255
256 // Send across data.
257 comm->Bcast(keys, 0);
258 comm->Bcast(vals, 0);
259 for (auto &cIt : m_compOrder)
260 {
261 comm->Bcast(cIt.second, 0);
262 }
263 m_meshGraph->SetCompositeOrdering(m_compOrder);
264
265 // Construct the keys and sizes of values for composite
266 // ordering
267 keys.resize(m_bndRegOrder.size());
268 vals.resize(m_bndRegOrder.size());
269
270 i = 0;
271 for (auto &bIt : m_bndRegOrder)
272 {
273 keys[i] = bIt.first;
274 vals[i++] = bIt.second.size();
275 }
276
277 // Send across data.
278 if (!keys.empty())
279 {
280 comm->Bcast(keys, 0);
281 }
282 if (!vals.empty())
283 {
284 comm->Bcast(vals, 0);
285 }
286 for (auto &bIt : m_bndRegOrder)
287 {
288 comm->Bcast(bIt.second, 0);
289 }
290 m_meshGraph->SetBndRegionOrdering(m_bndRegOrder);
291
292 if (session->DefinesCmdLineArgument("part-info"))
293 {
294 partitioner->PrintPartInfo(std::cout);
295 }
296 }
297 else
298 {
299 keys.resize(2);
300 comm->Bcast(keys, 0);
301
302 int cmpSize = keys[0];
303 int bndSize = keys[1];
304
305 keys.resize(cmpSize);
306 vals.resize(cmpSize);
307 comm->Bcast(keys, 0);
308 comm->Bcast(vals, 0);
309
310 for (int i = 0; i < keys.size(); ++i)
311 {
312 vector<unsigned int> tmp(vals[i]);
313 comm->Bcast(tmp, 0);
314 m_compOrder[keys[i]] = tmp;
315 }
316 m_meshGraph->SetCompositeOrdering(m_compOrder);
317
318 keys.resize(bndSize);
319 vals.resize(bndSize);
320 if (!keys.empty())
321 {
322 comm->Bcast(keys, 0);
323 }
324 if (!vals.empty())
325 {
326 comm->Bcast(vals, 0);
327 }
328 for (int i = 0; i < keys.size(); ++i)
329 {
330 vector<unsigned int> tmp(vals[i]);
331 comm->Bcast(tmp, 0);
332 m_bndRegOrder[keys[i]] = tmp;
333 }
334 m_meshGraph->SetBndRegionOrdering(m_bndRegOrder);
335 }
336 }
337 else
338 {
339 m_session->InitSession();
341
343 m_meshGraph->SetCompositeOrdering(m_compOrder);
344 auto comp = CreateCompositeDescriptor();
345
346 // Partitioner now operates in parallel. Each process receives
347 // partitioning over interconnect and writes its own session
348 // file to the working directory.
349 MeshPartitionSharedPtr partitioner =
351 partitionerName, session, session->GetComm(),
352 meshDimension, CreateMeshEntities(), comp);
353
354 partitioner->PartitionMesh(nParts, false);
355
356 vector<unsigned int> parts(1), tmp;
357 parts[0] = commMesh->GetRank();
358 vector<set<unsigned int>> elIDs(1);
359 partitioner->GetElementIDs(parts[0], tmp);
360 elIDs[0].insert(tmp.begin(), tmp.end());
361
362 // if (comm->GetTimeComm()->GetRank() == 0) // FIXME
363 // (OpenMPI 3.1.3)
364 {
365 this->WriteXMLGeometry(session->GetSessionName(), elIDs,
366 parts);
367 }
368
369 if (m_session->DefinesCmdLineArgument("part-info") && isRoot)
370 {
371 partitioner->PrintPartInfo(std::cout);
372 }
373 }
374
375 // Wait for all processors to finish their writing activities.
376 comm->Block();
377
378 std::string dirname = m_session->GetSessionName() + "_xml";
379 fs::path pdirname(dirname);
380 boost::format pad("P%1$07d.xml");
381 pad % comm->GetRowComm()->GetRank();
382 fs::path pFilename(pad.str());
383 fs::path fullpath = pdirname / pFilename;
384
385 std::vector<std::string> filenames = {
387 m_session->InitSession(filenames);
388 }
389 else if (!isRoot)
390 {
391 // No partitioning, non-root processors need to read the session
392 // file -- handles case where --npz is the same as number of
393 // processors.
394 m_session->InitSession();
395 }
396 }
397}
398
400 bool fillGraph)
401{
402 // Reset member variables.
403 m_meshGraph->Clear();
404
405 m_meshGraph->SetDomainRange(rng);
406 m_xmlGeom = m_session->GetElement("NEKTAR/GEOMETRY");
407
408 int err; /// Error value returned by TinyXML.
409
410 TiXmlAttribute *attr = m_xmlGeom->FirstAttribute();
411
412 // Initialize the mesh and space dimensions to 3 dimensions.
413 // We want to do this each time we read a file, so it should
414 // be done here and not just during class initialization.
415 int meshDimension = m_meshGraph->GetMeshDimension();
416 int spaceDimension = m_meshGraph->GetSpaceDimension();
417 m_meshPartitioned = false;
418 meshDimension = 3;
419 spaceDimension = 3;
420
421 while (attr)
422 {
423 std::string attrName(attr->Name());
424 if (attrName == "DIM")
425 {
426 err = attr->QueryIntValue(&meshDimension);
427 ASSERTL0(err == TIXML_SUCCESS, "Unable to read mesh dimension.");
428 }
429 else if (attrName == "SPACE")
430 {
431 err = attr->QueryIntValue(&spaceDimension);
432 ASSERTL0(err == TIXML_SUCCESS, "Unable to read space dimension.");
433 }
434 else if (attrName == "PARTITION")
435 {
436 err = attr->QueryIntValue(&m_partition);
437 ASSERTL0(err == TIXML_SUCCESS, "Unable to read partition.");
438 m_meshPartitioned = true;
439 m_meshGraph->SetMeshPartitioned(true);
440 }
441 else
442 {
443 std::string errstr("Unknown attribute: ");
444 errstr += attrName;
445 ASSERTL0(false, errstr.c_str());
446 }
447
448 // Get the next attribute.
449 attr = attr->Next();
450 }
451
452 m_meshGraph->SetMeshDimension(meshDimension);
453 m_meshGraph->SetSpaceDimension(spaceDimension);
454
455 ASSERTL0(meshDimension <= spaceDimension,
456 "Mesh dimension greater than space dimension");
457
459 v_ReadCurves();
460 if (meshDimension >= 2)
461 {
462 v_ReadEdges();
463 if (meshDimension == 3)
464 {
465 v_ReadFaces();
466 }
467 }
468 ReadElements();
470 ReadDomain();
471
472 if (fillGraph)
473 {
474 m_meshGraph->FillGraph();
475 }
476}
477
479{
480 auto &vertSet = m_meshGraph->GetAllPointGeoms();
481 int spaceDimension = m_meshGraph->GetSpaceDimension();
482
483 // Now read the vertices
484 TiXmlElement *element = m_xmlGeom->FirstChildElement("VERTEX");
485 ASSERTL0(element, "Unable to find mesh VERTEX tag in file.");
486
487 NekDouble xscale, yscale, zscale;
488
489 // check to see if any scaling parameters are in
490 // attributes and determine these values
491 LibUtilities::Interpreter expEvaluator;
492 const char *xscal = element->Attribute("XSCALE");
493 if (!xscal)
494 {
495 xscale = 1.0;
496 }
497 else
498 {
499 std::string xscalstr = xscal;
500 int expr_id = expEvaluator.DefineFunction("", xscalstr);
501 xscale = expEvaluator.Evaluate(expr_id);
502 }
503
504 const char *yscal = element->Attribute("YSCALE");
505 if (!yscal)
506 {
507 yscale = 1.0;
508 }
509 else
510 {
511 std::string yscalstr = yscal;
512 int expr_id = expEvaluator.DefineFunction("", yscalstr);
513 yscale = expEvaluator.Evaluate(expr_id);
514 }
515
516 const char *zscal = element->Attribute("ZSCALE");
517 if (!zscal)
518 {
519 zscale = 1.0;
520 }
521 else
522 {
523 std::string zscalstr = zscal;
524 int expr_id = expEvaluator.DefineFunction("", zscalstr);
525 zscale = expEvaluator.Evaluate(expr_id);
526 }
527
528 NekDouble xmove, ymove, zmove;
529
530 // check to see if any moving parameters are in
531 // attributes and determine these values
532
533 const char *xmov = element->Attribute("XMOVE");
534 if (!xmov)
535 {
536 xmove = 0.0;
537 }
538 else
539 {
540 std::string xmovstr = xmov;
541 int expr_id = expEvaluator.DefineFunction("", xmovstr);
542 xmove = expEvaluator.Evaluate(expr_id);
543 }
544
545 const char *ymov = element->Attribute("YMOVE");
546 if (!ymov)
547 {
548 ymove = 0.0;
549 }
550 else
551 {
552 std::string ymovstr = ymov;
553 int expr_id = expEvaluator.DefineFunction("", ymovstr);
554 ymove = expEvaluator.Evaluate(expr_id);
555 }
556
557 const char *zmov = element->Attribute("ZMOVE");
558 if (!zmov)
559 {
560 zmove = 0.0;
561 }
562 else
563 {
564 std::string zmovstr = zmov;
565 int expr_id = expEvaluator.DefineFunction("", zmovstr);
566 zmove = expEvaluator.Evaluate(expr_id);
567 }
568
569 TiXmlElement *vertex = element->FirstChildElement("V");
570
571 int indx;
572
573 while (vertex)
574 {
575 TiXmlAttribute *vertexAttr = vertex->FirstAttribute();
576 std::string attrName(vertexAttr->Name());
577
578 ASSERTL0(attrName == "ID",
579 (std::string("Unknown attribute name: ") + attrName).c_str());
580
581 int err = vertexAttr->QueryIntValue(&indx);
582 ASSERTL0(err == TIXML_SUCCESS, "Unable to read attribute ID.");
583
584 // Now read body of vertex
585 std::string vertexBodyStr;
586
587 TiXmlNode *vertexBody = vertex->FirstChild();
588
589 while (vertexBody)
590 {
591 // Accumulate all non-comment body data.
592 if (vertexBody->Type() == TiXmlNode::TINYXML_TEXT)
593 {
594 vertexBodyStr += vertexBody->ToText()->Value();
595 vertexBodyStr += " ";
596 }
597
598 vertexBody = vertexBody->NextSibling();
599 }
600
601 ASSERTL0(!vertexBodyStr.empty(),
602 "Vertex definitions must contain vertex data.");
603
604 // Get vertex data from the data string.
605 NekDouble xval, yval, zval;
606 std::istringstream vertexDataStrm(vertexBodyStr.c_str());
607
608 try
609 {
610 while (!vertexDataStrm.fail())
611 {
612 vertexDataStrm >> xval >> yval >> zval;
613
614 xval = xval * xscale + xmove;
615 yval = yval * yscale + ymove;
616 zval = zval * zscale + zmove;
617
618 // Need to check it here because we may not be
619 // good after the read indicating that there
620 // was nothing to read.
621 if (!vertexDataStrm.fail())
622 {
625 spaceDimension, indx, xval, yval, zval));
626 vertSet[indx] = vert;
627 }
628 }
629 }
630 catch (...)
631 {
632 ASSERTL0(false, "Unable to read VERTEX data.");
633 }
634
635 vertex = vertex->NextSiblingElement("V");
636 }
637}
638
640{
641 auto &curvedEdges = m_meshGraph->GetCurvedEdges();
642 auto &curvedFaces = m_meshGraph->GetCurvedFaces();
643 int meshDimension = m_meshGraph->GetMeshDimension();
644
645 // check to see if any scaling parameters are in
646 // attributes and determine these values
647 TiXmlElement *element = m_xmlGeom->FirstChildElement("VERTEX");
648 ASSERTL0(element, "Unable to find mesh VERTEX tag in file.");
649
650 NekDouble xscale, yscale, zscale;
651
652 LibUtilities::Interpreter expEvaluator;
653 const char *xscal = element->Attribute("XSCALE");
654 if (!xscal)
655 {
656 xscale = 1.0;
657 }
658 else
659 {
660 std::string xscalstr = xscal;
661 int expr_id = expEvaluator.DefineFunction("", xscalstr);
662 xscale = expEvaluator.Evaluate(expr_id);
663 }
664
665 const char *yscal = element->Attribute("YSCALE");
666 if (!yscal)
667 {
668 yscale = 1.0;
669 }
670 else
671 {
672 std::string yscalstr = yscal;
673 int expr_id = expEvaluator.DefineFunction("", yscalstr);
674 yscale = expEvaluator.Evaluate(expr_id);
675 }
676
677 const char *zscal = element->Attribute("ZSCALE");
678 if (!zscal)
679 {
680 zscale = 1.0;
681 }
682 else
683 {
684 std::string zscalstr = zscal;
685 int expr_id = expEvaluator.DefineFunction("", zscalstr);
686 zscale = expEvaluator.Evaluate(expr_id);
687 }
688
689 NekDouble xmove, ymove, zmove;
690
691 // check to see if any moving parameters are in
692 // attributes and determine these values
693
694 const char *xmov = element->Attribute("XMOVE");
695 if (!xmov)
696 {
697 xmove = 0.0;
698 }
699 else
700 {
701 std::string xmovstr = xmov;
702 int expr_id = expEvaluator.DefineFunction("", xmovstr);
703 xmove = expEvaluator.Evaluate(expr_id);
704 }
705
706 const char *ymov = element->Attribute("YMOVE");
707 if (!ymov)
708 {
709 ymove = 0.0;
710 }
711 else
712 {
713 std::string ymovstr = ymov;
714 int expr_id = expEvaluator.DefineFunction("", ymovstr);
715 ymove = expEvaluator.Evaluate(expr_id);
716 }
717
718 const char *zmov = element->Attribute("ZMOVE");
719 if (!zmov)
720 {
721 zmove = 0.0;
722 }
723 else
724 {
725 std::string zmovstr = zmov;
726 int expr_id = expEvaluator.DefineFunction("", zmovstr);
727 zmove = expEvaluator.Evaluate(expr_id);
728 }
729
730 int err;
731
732 /// Look for elements in CURVE block.
733 TiXmlElement *field = m_xmlGeom->FirstChildElement("CURVED");
734
735 if (!field) // return if no curved entities
736 {
737 return;
738 }
739
740 /// All curves are of the form: "<? ID="#" TYPE="GLL OR other
741 /// points type" NUMPOINTS="#"> ... </?>", with ? being an
742 /// element type (either E or F).
743
744 TiXmlElement *edgelement = field->FirstChildElement("E");
745
746 int edgeindx, edgeid;
747
748 while (edgelement)
749 {
750 std::string edge(edgelement->ValueStr());
751 ASSERTL0(edge == "E",
752 (std::string("Unknown 3D curve type:") + edge).c_str());
753
754 /// Read id attribute.
755 err = edgelement->QueryIntAttribute("ID", &edgeindx);
756 ASSERTL0(err == TIXML_SUCCESS, "Unable to read curve attribute ID.");
757
758 /// Read edge id attribute.
759 err = edgelement->QueryIntAttribute("EDGEID", &edgeid);
760 ASSERTL0(err == TIXML_SUCCESS,
761 "Unable to read curve attribute EDGEID.");
762
763 /// Read text edgelement description.
764 std::string elementStr;
765 TiXmlNode *elementChild = edgelement->FirstChild();
766
767 while (elementChild)
768 {
769 // Accumulate all non-comment element data
770 if (elementChild->Type() == TiXmlNode::TINYXML_TEXT)
771 {
772 elementStr += elementChild->ToText()->ValueStr();
773 elementStr += " ";
774 }
775 elementChild = elementChild->NextSibling();
776 }
777
778 ASSERTL0(!elementStr.empty(), "Unable to read curve description body.");
779
780 /// Parse out the element components corresponding to type of
781 /// element.
782 if (edge == "E")
783 {
784 int numPts = 0;
785 // Determine the points type
786 std::string typeStr = edgelement->Attribute("TYPE");
787 ASSERTL0(!typeStr.empty(), "TYPE must be specified in "
788 "points definition");
789
791 const std::string *begStr = LibUtilities::kPointsTypeStr;
792 const std::string *endStr =
794 const std::string *ptsStr = std::find(begStr, endStr, typeStr);
795
796 ASSERTL0(ptsStr != endStr, "Invalid points type.");
797 type = (LibUtilities::PointsType)(ptsStr - begStr);
798
799 // Determine the number of points
800 err = edgelement->QueryIntAttribute("NUMPOINTS", &numPts);
801 ASSERTL0(err == TIXML_SUCCESS,
802 "Unable to read curve attribute NUMPOINTS.");
803 CurveSharedPtr curve(
805
806 // Read points (x, y, z)
807 NekDouble xval, yval, zval;
808 std::istringstream elementDataStrm(elementStr.c_str());
809 try
810 {
811 while (!elementDataStrm.fail())
812 {
813 elementDataStrm >> xval >> yval >> zval;
814
815 xval = xval * xscale + xmove;
816 yval = yval * yscale + ymove;
817 zval = zval * zscale + zmove;
818
819 // Need to check it here because we may not be
820 // good after the read indicating that there
821 // was nothing to read.
822 if (!elementDataStrm.fail())
823 {
826 meshDimension, edgeindx, xval, yval, zval));
827
828 curve->m_points.push_back(vert);
829 }
830 }
831 }
832 catch (...)
833 {
835 (std::string("Unable to read curve data for EDGE: ") +
836 elementStr)
837 .c_str());
838 }
839
840 ASSERTL0(curve->m_points.size() == numPts,
841 "Number of points specificed by attribute "
842 "NUMPOINTS is different from number of points "
843 "in list (edgeid = " +
844 boost::lexical_cast<string>(edgeid));
845
846 curvedEdges[edgeid] = curve;
847
848 edgelement = edgelement->NextSiblingElement("E");
849
850 } // end if-loop
851
852 } // end while-loop
853
854 TiXmlElement *facelement = field->FirstChildElement("F");
855 int faceindx, faceid;
856
857 while (facelement)
858 {
859 std::string face(facelement->ValueStr());
860 ASSERTL0(face == "F",
861 (std::string("Unknown 3D curve type: ") + face).c_str());
862
863 /// Read id attribute.
864 err = facelement->QueryIntAttribute("ID", &faceindx);
865 ASSERTL0(err == TIXML_SUCCESS, "Unable to read curve attribute ID.");
866
867 /// Read face id attribute.
868 err = facelement->QueryIntAttribute("FACEID", &faceid);
869 ASSERTL0(err == TIXML_SUCCESS,
870 "Unable to read curve attribute FACEID.");
871
872 /// Read text face element description.
873 std::string elementStr;
874 TiXmlNode *elementChild = facelement->FirstChild();
875
876 while (elementChild)
877 {
878 // Accumulate all non-comment element data
879 if (elementChild->Type() == TiXmlNode::TINYXML_TEXT)
880 {
881 elementStr += elementChild->ToText()->ValueStr();
882 elementStr += " ";
883 }
884 elementChild = elementChild->NextSibling();
885 }
886
887 ASSERTL0(!elementStr.empty(), "Unable to read curve description body.");
888
889 /// Parse out the element components corresponding to type of
890 /// element.
891 if (face == "F")
892 {
893 std::string typeStr = facelement->Attribute("TYPE");
894 ASSERTL0(!typeStr.empty(), "TYPE must be specified in "
895 "points definition");
897 const std::string *begStr = LibUtilities::kPointsTypeStr;
898 const std::string *endStr =
900 const std::string *ptsStr = std::find(begStr, endStr, typeStr);
901
902 ASSERTL0(ptsStr != endStr, "Invalid points type.");
903 type = (LibUtilities::PointsType)(ptsStr - begStr);
904
905 std::string numptsStr = facelement->Attribute("NUMPOINTS");
906 ASSERTL0(!numptsStr.empty(),
907 "NUMPOINTS must be specified in points definition");
908 int numPts = 0;
909 std::stringstream s;
910 s << numptsStr;
911 s >> numPts;
912
913 CurveSharedPtr curve(
915
916 ASSERTL0(numPts >= 3, "NUMPOINTS for face must be greater than 2");
917
918 if (numPts == 3)
919 {
920 ASSERTL0(ptsStr != endStr, "Invalid points type.");
921 }
922
923 // Read points (x, y, z)
924 NekDouble xval, yval, zval;
925 std::istringstream elementDataStrm(elementStr.c_str());
926 try
927 {
928 while (!elementDataStrm.fail())
929 {
930 elementDataStrm >> xval >> yval >> zval;
931
932 // Need to check it here because we
933 // may not be good after the read
934 // indicating that there was nothing
935 // to read.
936 if (!elementDataStrm.fail())
937 {
940 meshDimension, faceindx, xval, yval, zval));
941 curve->m_points.push_back(vert);
942 }
943 }
944 }
945 catch (...)
946 {
948 (std::string("Unable to read curve data for FACE: ") +
949 elementStr)
950 .c_str());
951 }
952 curvedFaces[faceid] = curve;
953
954 facelement = facelement->NextSiblingElement("F");
955 }
956 }
957}
958
960{
961 auto &domain = m_meshGraph->GetDomain();
962
963 TiXmlElement *element = nullptr;
964 /// Look for data in DOMAIN block.
965 element = m_xmlGeom->FirstChildElement("DOMAIN");
966
967 ASSERTL0(element, "Unable to find DOMAIN tag in file.");
968
969 /// Elements are of the form: "<D ID = "N"> ... </D>".
970 /// Read the ID field first.
971 TiXmlElement *multidomains = element->FirstChildElement("D");
972
973 if (multidomains)
974 {
975 while (multidomains)
976 {
977 int indx;
978 int err = multidomains->QueryIntAttribute("ID", &indx);
979 ASSERTL0(err == TIXML_SUCCESS,
980 "Unable to read attribute ID in Domain.");
981
982 TiXmlNode *elementChild = multidomains->FirstChild();
983 while (elementChild &&
984 elementChild->Type() != TiXmlNode::TINYXML_TEXT)
985 {
986 elementChild = elementChild->NextSibling();
987 }
988
989 ASSERTL0(elementChild, "Unable to read DOMAIN body.");
990 std::string elementStr = elementChild->ToText()->ValueStr();
991
992 elementStr = elementStr.substr(elementStr.find_first_not_of(" "));
993
994 std::string::size_type indxBeg = elementStr.find_first_of('[') + 1;
995 std::string::size_type indxEnd = elementStr.find_last_of(']') - 1;
996 std::string indxStr =
997 elementStr.substr(indxBeg, indxEnd - indxBeg + 1);
998
999 ASSERTL0(
1000 !indxStr.empty(),
1001 "Unable to read domain's composite index (index missing?).");
1002
1003 // Read the domain composites.
1004 // Parse the composites into a list.
1005 map<int, CompositeSharedPtr> unrollDomain;
1006 m_meshGraph->GetCompositeList(indxStr, unrollDomain);
1007 domain[indx] = unrollDomain;
1008
1009 ASSERTL0(!domain[indx].empty(),
1010 (std::string(
1011 "Unable to obtain domain's referenced composite: ") +
1012 indxStr)
1013 .c_str());
1014
1015 /// Keep looking
1016 multidomains = multidomains->NextSiblingElement("D");
1017 }
1018 }
1019 else // previous definition of just one composite
1020 {
1021
1022 // find the non comment portion of the body.
1023 TiXmlNode *elementChild = element->FirstChild();
1024 while (elementChild && elementChild->Type() != TiXmlNode::TINYXML_TEXT)
1025 {
1026 elementChild = elementChild->NextSibling();
1027 }
1028
1029 ASSERTL0(elementChild, "Unable to read DOMAIN body.");
1030 std::string elementStr = elementChild->ToText()->ValueStr();
1031
1032 elementStr = elementStr.substr(elementStr.find_first_not_of(" "));
1033
1034 std::string::size_type indxBeg = elementStr.find_first_of('[') + 1;
1035 std::string::size_type indxEnd = elementStr.find_last_of(']') - 1;
1036 std::string indxStr = elementStr.substr(indxBeg, indxEnd - indxBeg + 1);
1037
1038 ASSERTL0(!indxStr.empty(),
1039 "Unable to read domain's composite index (index missing?).");
1040
1041 // Read the domain composites.
1042 // Parse the composites into a list.
1043 map<int, CompositeSharedPtr> fullDomain;
1044 m_meshGraph->GetCompositeList(indxStr, fullDomain);
1045 domain[0] = fullDomain;
1046
1047 ASSERTL0(
1048 !domain[0].empty(),
1049 (std::string("Unable to obtain domain's referenced composite: ") +
1050 indxStr)
1051 .c_str());
1052 }
1053}
1054
1056{
1057 auto &segGeoms = m_meshGraph->GetAllSegGeoms();
1058 auto &curvedEdges = m_meshGraph->GetCurvedEdges();
1059 int spaceDimension = m_meshGraph->GetSpaceDimension();
1060
1061 CurveMap::iterator it;
1062
1063 /// Look for elements in ELEMENT block.
1064 TiXmlElement *field = m_xmlGeom->FirstChildElement("EDGE");
1065
1066 ASSERTL0(field, "Unable to find EDGE tag in file.");
1067
1068 /// All elements are of the form: "<E ID="#"> ... </E>", with
1069 /// ? being the element type.
1070 /// Read the ID field first.
1071 TiXmlElement *edgeElement = field->FirstChildElement("E");
1072
1073 /// Since all edge data is one big text block, we need to accumulate
1074 /// all TINYXML_TEXT data and then parse it. This approach effectively
1075 /// skips
1076 /// all comments or other node types since we only care about the
1077 /// edge list. We cannot handle missing edge numbers as we could
1078 /// with missing element numbers due to the text block format.
1079 std::string edgeStr;
1080 int indx;
1081
1082 while (edgeElement)
1083 {
1084 int err = edgeElement->QueryIntAttribute("ID", &indx);
1085 ASSERTL0(err == TIXML_SUCCESS, "Unable to read edge attribute ID.");
1086
1087 TiXmlNode *child = edgeElement->FirstChild();
1088 edgeStr.clear();
1089 if (child->Type() == TiXmlNode::TINYXML_TEXT)
1090 {
1091 edgeStr += child->ToText()->ValueStr();
1092 }
1093
1094 /// Now parse out the edges, three fields at a time.
1095 int vertex1, vertex2;
1096 std::istringstream edgeDataStrm(edgeStr.c_str());
1097
1098 try
1099 {
1100 while (!edgeDataStrm.fail())
1101 {
1102 edgeDataStrm >> vertex1 >> vertex2;
1103
1104 // Must check after the read because we
1105 // may be at the end and not know it. If
1106 // we are at the end we will add a
1107 // duplicate of the last entry if we don't
1108 // check here.
1109 if (!edgeDataStrm.fail())
1110 {
1111 PointGeomSharedPtr vertices[2] = {
1112 m_meshGraph->GetVertex(vertex1),
1113 m_meshGraph->GetVertex(vertex2)};
1114 it = curvedEdges.find(indx);
1115
1116 if (it == curvedEdges.end())
1117 {
1118 segGeoms[indx] =
1120 indx, spaceDimension, vertices);
1121 }
1122 else
1123 {
1124 segGeoms[indx] =
1126 indx, spaceDimension, vertices, it->second);
1127 }
1128 }
1129 }
1130 }
1131 catch (...)
1132 {
1133 NEKERROR(
1135 (std::string("Unable to read edge data: ") + edgeStr).c_str());
1136 }
1137
1138 edgeElement = edgeElement->NextSiblingElement("E");
1139 }
1140}
1141
1143{
1144 auto &curvedFaces = m_meshGraph->GetCurvedFaces();
1145 auto &triGeoms = m_meshGraph->GetAllTriGeoms();
1146 auto &quadGeoms = m_meshGraph->GetAllQuadGeoms();
1147
1148 /// Look for elements in FACE block.
1149 TiXmlElement *field = m_xmlGeom->FirstChildElement("FACE");
1150
1151 ASSERTL0(field, "Unable to find FACE tag in file.");
1152
1153 /// All faces are of the form: "<? ID="#"> ... </?>", with
1154 /// ? being an element type (either Q or T).
1155 /// They might be in compressed format and so then need upacking.
1156
1157 TiXmlElement *element = field->FirstChildElement();
1158 CurveMap::iterator it;
1159
1160 while (element)
1161 {
1162 std::string elementType(element->ValueStr());
1163
1164 ASSERTL0(elementType == "Q" || elementType == "T",
1165 (std::string("Unknown 3D face type: ") + elementType).c_str());
1166
1167 /// Read id attribute.
1168 int indx;
1169 int err = element->QueryIntAttribute("ID", &indx);
1170 ASSERTL0(err == TIXML_SUCCESS, "Unable to read face attribute ID.");
1171
1172 /// See if this face has curves.
1173 it = curvedFaces.find(indx);
1174
1175 /// Read text element description.
1176 TiXmlNode *elementChild = element->FirstChild();
1177 std::string elementStr;
1178 while (elementChild)
1179 {
1180 if (elementChild->Type() == TiXmlNode::TINYXML_TEXT)
1181 {
1182 elementStr += elementChild->ToText()->ValueStr();
1183 }
1184 elementChild = elementChild->NextSibling();
1185 }
1186
1187 ASSERTL0(!elementStr.empty(), "Unable to read face description body.");
1188
1189 /// Parse out the element components corresponding to type of
1190 /// element.
1191 if (elementType == "T")
1192 {
1193 // Read three edge numbers
1194 int edge1, edge2, edge3;
1195 std::istringstream elementDataStrm(elementStr.c_str());
1196
1197 try
1198 {
1199 elementDataStrm >> edge1;
1200 elementDataStrm >> edge2;
1201 elementDataStrm >> edge3;
1202
1203 ASSERTL0(
1204 !elementDataStrm.fail(),
1205 (std::string("Unable to read face data for TRIANGLE: ") +
1206 elementStr)
1207 .c_str());
1208
1209 /// Create a TriGeom to hold the new definition.
1211 m_meshGraph->GetSegGeom(edge1),
1212 m_meshGraph->GetSegGeom(edge2),
1213 m_meshGraph->GetSegGeom(edge3)};
1214
1215 if (it == curvedFaces.end())
1216 {
1217 triGeoms[indx] =
1219 }
1220 else
1221 {
1223 indx, edges, it->second);
1224 }
1225 triGeoms[indx]->SetGlobalID(indx);
1226 }
1227 catch (...)
1228 {
1229 NEKERROR(
1231 (std::string("Unable to read face data for TRIANGLE: ") +
1232 elementStr)
1233 .c_str());
1234 }
1235 }
1236 else if (elementType == "Q")
1237 {
1238 // Read four edge numbers
1239 int edge1, edge2, edge3, edge4;
1240 std::istringstream elementDataStrm(elementStr.c_str());
1241
1242 try
1243 {
1244 elementDataStrm >> edge1;
1245 elementDataStrm >> edge2;
1246 elementDataStrm >> edge3;
1247 elementDataStrm >> edge4;
1248
1249 ASSERTL0(!elementDataStrm.fail(),
1250 (std::string("Unable to read face data for QUAD: ") +
1251 elementStr)
1252 .c_str());
1253
1254 /// Create a QuadGeom to hold the new definition.
1256 m_meshGraph->GetSegGeom(edge1),
1257 m_meshGraph->GetSegGeom(edge2),
1258 m_meshGraph->GetSegGeom(edge3),
1259 m_meshGraph->GetSegGeom(edge4)};
1260
1261 QuadGeomSharedPtr quadgeom;
1262
1263 if (it == curvedFaces.end())
1264 {
1265 quadGeoms[indx] =
1267 }
1268 else
1269 {
1270 quadGeoms[indx] =
1272 it->second);
1273 }
1274 quadGeoms[indx]->SetGlobalID(indx);
1275 }
1276 catch (...)
1277 {
1279 (std::string("Unable to read face data for QUAD: ") +
1280 elementStr)
1281 .c_str());
1282 }
1283 }
1284 // Keep looking
1285 element = element->NextSiblingElement();
1286 }
1287}
1288
1290{
1291 int meshDimension = m_meshGraph->GetMeshDimension();
1292
1293 switch (meshDimension)
1294 {
1295 case 1:
1297 break;
1298 case 2:
1300 break;
1301 case 3:
1303 break;
1304 }
1305}
1306
1308{
1309 auto &curvedEdges = m_meshGraph->GetCurvedEdges();
1310 auto &segGeoms = m_meshGraph->GetAllSegGeoms();
1311 int spaceDimension = m_meshGraph->GetSpaceDimension();
1312
1313 TiXmlElement *field = nullptr;
1314
1315 /// Look for elements in ELEMENT block.
1316 field = m_xmlGeom->FirstChildElement("ELEMENT");
1317
1318 ASSERTL0(field, "Unable to find ELEMENT tag in file.");
1319
1320 /// All elements are of the form: "<S ID = n> ... </S>", with
1321 /// ? being the element type.
1322
1323 TiXmlElement *segment = field->FirstChildElement("S");
1324 CurveMap::iterator it;
1325
1326 while (segment)
1327 {
1328 int indx;
1329 int err = segment->QueryIntAttribute("ID", &indx);
1330 ASSERTL0(err == TIXML_SUCCESS, "Unable to read element attribute ID.");
1331
1332 TiXmlNode *elementChild = segment->FirstChild();
1333 while (elementChild && elementChild->Type() != TiXmlNode::TINYXML_TEXT)
1334 {
1335 elementChild = elementChild->NextSibling();
1336 }
1337
1338 ASSERTL0(elementChild, "Unable to read element description body.");
1339 std::string elementStr = elementChild->ToText()->ValueStr();
1340
1341 /// Parse out the element components corresponding to type of
1342 /// element.
1343 /// Read two vertex numbers
1344 int vertex1, vertex2;
1345 std::istringstream elementDataStrm(elementStr.c_str());
1346
1347 try
1348 {
1349 elementDataStrm >> vertex1;
1350 elementDataStrm >> vertex2;
1351
1352 ASSERTL0(!elementDataStrm.fail(),
1353 (std::string("Unable to read element data for SEGMENT: ") +
1354 elementStr)
1355 .c_str());
1356
1357 PointGeomSharedPtr vertices[2] = {m_meshGraph->GetVertex(vertex1),
1358 m_meshGraph->GetVertex(vertex2)};
1359 it = curvedEdges.find(indx);
1360
1361 if (it == curvedEdges.end())
1362 {
1364 indx, spaceDimension, vertices);
1365 }
1366 else
1367 {
1369 indx, spaceDimension, vertices, it->second);
1370 }
1371 segGeoms[indx]->SetGlobalID(indx);
1372 }
1373 catch (...)
1374 {
1376 (std::string("Unable to read element data for segment: ") +
1377 elementStr)
1378 .c_str());
1379 }
1380 /// Keep looking for additional segments
1381 segment = segment->NextSiblingElement("S");
1382 }
1383}
1384
1386{
1387 auto &curvedFaces = m_meshGraph->GetCurvedFaces();
1388 auto &triGeoms = m_meshGraph->GetAllTriGeoms();
1389 auto &quadGeoms = m_meshGraph->GetAllQuadGeoms();
1390
1391 /// Look for elements in ELEMENT block.
1392 TiXmlElement *field = m_xmlGeom->FirstChildElement("ELEMENT");
1393
1394 ASSERTL0(field, "Unable to find ELEMENT tag in file.");
1395
1396 // Set up curve map for curved elements on an embedded manifold.
1397 CurveMap::iterator it;
1398
1399 /// All elements are of the form: "<? ID="#"> ... </?>", with
1400 /// ? being the element type.
1401
1402 TiXmlElement *element = field->FirstChildElement();
1403
1404 while (element)
1405 {
1406 std::string elementType(element->ValueStr());
1407
1408 ASSERTL0(
1409 elementType == "Q" || elementType == "T",
1410 (std::string("Unknown 2D element type: ") + elementType).c_str());
1411
1412 /// Read id attribute.
1413 int indx;
1414 int err = element->QueryIntAttribute("ID", &indx);
1415 ASSERTL0(err == TIXML_SUCCESS, "Unable to read element attribute ID.");
1416
1417 it = curvedFaces.find(indx);
1418
1419 /// Read text element description.
1420 TiXmlNode *elementChild = element->FirstChild();
1421 std::string elementStr;
1422 while (elementChild)
1423 {
1424 if (elementChild->Type() == TiXmlNode::TINYXML_TEXT)
1425 {
1426 elementStr += elementChild->ToText()->ValueStr();
1427 }
1428 elementChild = elementChild->NextSibling();
1429 }
1430
1431 ASSERTL0(!elementStr.empty(),
1432 "Unable to read element description body.");
1433
1434 /// Parse out the element components corresponding to type of
1435 /// element.
1436 if (elementType == "T")
1437 {
1438 // Read three edge numbers
1439 int edge1, edge2, edge3;
1440 std::istringstream elementDataStrm(elementStr.c_str());
1441
1442 try
1443 {
1444 elementDataStrm >> edge1;
1445 elementDataStrm >> edge2;
1446 elementDataStrm >> edge3;
1447
1448 ASSERTL0(
1449 !elementDataStrm.fail(),
1450 (std::string("Unable to read element data for TRIANGLE: ") +
1451 elementStr)
1452 .c_str());
1453
1454 /// Create a TriGeom to hold the new definition.
1456 m_meshGraph->GetSegGeom(edge1),
1457 m_meshGraph->GetSegGeom(edge2),
1458 m_meshGraph->GetSegGeom(edge3)};
1459
1460 if (it == curvedFaces.end())
1461 {
1462 triGeoms[indx] =
1464 }
1465 else
1466 {
1468 indx, edges, it->second);
1469 }
1470 triGeoms[indx]->SetGlobalID(indx);
1471 }
1472 catch (...)
1473 {
1474 NEKERROR(
1476 (std::string("Unable to read element data for TRIANGLE: ") +
1477 elementStr)
1478 .c_str());
1479 }
1480 }
1481 else if (elementType == "Q")
1482 {
1483 // Read four edge numbers
1484 int edge1, edge2, edge3, edge4;
1485 std::istringstream elementDataStrm(elementStr.c_str());
1486
1487 try
1488 {
1489 elementDataStrm >> edge1;
1490 elementDataStrm >> edge2;
1491 elementDataStrm >> edge3;
1492 elementDataStrm >> edge4;
1493
1494 ASSERTL0(
1495 !elementDataStrm.fail(),
1496 (std::string("Unable to read element data for QUAD: ") +
1497 elementStr)
1498 .c_str());
1499
1500 /// Create a QuadGeom to hold the new definition.
1502 m_meshGraph->GetSegGeom(edge1),
1503 m_meshGraph->GetSegGeom(edge2),
1504 m_meshGraph->GetSegGeom(edge3),
1505 m_meshGraph->GetSegGeom(edge4)};
1506
1507 QuadGeomSharedPtr quadgeom;
1508 if (it == curvedFaces.end())
1509 {
1510 quadGeoms[indx] =
1512 }
1513 else
1514 {
1515 quadGeoms[indx] =
1517 it->second);
1518 }
1519 quadGeoms[indx]->SetGlobalID(indx);
1520 }
1521 catch (...)
1522 {
1523 NEKERROR(
1525 (std::string("Unable to read element data for QUAD: ") +
1526 elementStr)
1527 .c_str());
1528 }
1529 }
1530 /// Keep looking
1531 element = element->NextSiblingElement();
1532 }
1533}
1534
1536{
1537 auto &tetGeoms = m_meshGraph->GetAllTetGeoms();
1538 auto &pyrGeoms = m_meshGraph->GetAllPyrGeoms();
1539 auto &prismGeoms = m_meshGraph->GetAllPrismGeoms();
1540 auto &hexGeoms = m_meshGraph->GetAllHexGeoms();
1541
1542 /// Look for elements in ELEMENT block.
1543 TiXmlElement *field = m_xmlGeom->FirstChildElement("ELEMENT");
1544
1545 ASSERTL0(field, "Unable to find ELEMENT tag in file.");
1546
1547 /// All elements are of the form: "<? ID="#"> ... </?>", with
1548 /// ? being the element type.
1549
1550 TiXmlElement *element = field->FirstChildElement();
1551
1552 while (element)
1553 {
1554 std::string elementType(element->ValueStr());
1555
1556 // A - tet, P - pyramid, R - prism, H - hex
1557 ASSERTL0(
1558 elementType == "A" || elementType == "P" || elementType == "R" ||
1559 elementType == "H",
1560 (std::string("Unknown 3D element type: ") + elementType).c_str());
1561
1562 /// Read id attribute.
1563 int indx;
1564 int err = element->QueryIntAttribute("ID", &indx);
1565 ASSERTL0(err == TIXML_SUCCESS, "Unable to read element attribute ID.");
1566
1567 /// Read text element description.
1568 TiXmlNode *elementChild = element->FirstChild();
1569 std::string elementStr;
1570 while (elementChild)
1571 {
1572 if (elementChild->Type() == TiXmlNode::TINYXML_TEXT)
1573 {
1574 elementStr += elementChild->ToText()->ValueStr();
1575 }
1576 elementChild = elementChild->NextSibling();
1577 }
1578
1579 ASSERTL0(!elementStr.empty(),
1580 "Unable to read element description body.");
1581
1582 std::istringstream elementDataStrm(elementStr.c_str());
1583
1584 /// Parse out the element components corresponding to type of
1585 /// element.
1586
1587 // Tetrahedral
1588 if (elementType == "A")
1589 {
1590 try
1591 {
1592 /// Create arrays for the tri and quad faces.
1593 const int kNfaces = TetGeom::kNfaces;
1594 const int kNtfaces = TetGeom::kNtfaces;
1595 const int kNqfaces = TetGeom::kNqfaces;
1596 TriGeomSharedPtr tfaces[kNtfaces];
1597 int Ntfaces = 0;
1598 int Nqfaces = 0;
1599
1600 /// Fill the arrays and make sure there aren't too many
1601 /// faces.
1602 std::stringstream errorstring;
1603 errorstring << "Element " << indx << " must have " << kNtfaces
1604 << " triangle face(s), and " << kNqfaces
1605 << " quadrilateral face(s).";
1606 for (int i = 0; i < kNfaces; i++)
1607 {
1608 int faceID;
1609 elementDataStrm >> faceID;
1610 Geometry2DSharedPtr face =
1611 m_meshGraph->GetGeometry2D(faceID);
1612 if (face == Geometry2DSharedPtr() ||
1613 (face->GetShapeType() != LibUtilities::eTriangle &&
1614 face->GetShapeType() != LibUtilities::eQuadrilateral))
1615 {
1616 std::stringstream errorstring;
1617 errorstring << "Element " << indx
1618 << " has invalid face: " << faceID;
1619 ASSERTL0(false, errorstring.str().c_str());
1620 }
1621 else if (face->GetShapeType() == LibUtilities::eTriangle)
1622 {
1623 ASSERTL0(Ntfaces < kNtfaces, errorstring.str().c_str());
1624 tfaces[Ntfaces++] = static_pointer_cast<TriGeom>(face);
1625 }
1626 else if (face->GetShapeType() ==
1628 {
1629 ASSERTL0(Nqfaces < kNqfaces, errorstring.str().c_str());
1630 }
1631 }
1632
1633 /// Make sure all of the face indicies could be read, and
1634 /// that there weren't too few.
1635 ASSERTL0(!elementDataStrm.fail(),
1636 (std::string(
1637 "Unable to read element data for TETRAHEDRON: ") +
1638 elementStr)
1639 .c_str());
1640 ASSERTL0(Ntfaces == kNtfaces, errorstring.str().c_str());
1641 ASSERTL0(Nqfaces == kNqfaces, errorstring.str().c_str());
1642
1643 TetGeomSharedPtr tetgeom(
1645
1646 tetGeoms[indx] = tetgeom;
1647 m_meshGraph->PopulateFaceToElMap(tetgeom, kNfaces);
1648 }
1649 catch (...)
1650 {
1652 (std::string(
1653 "Unable to read element data for TETRAHEDRON: ") +
1654 elementStr)
1655 .c_str());
1656 }
1657 }
1658 // Pyramid
1659 else if (elementType == "P")
1660 {
1661 try
1662 {
1663 /// Create arrays for the tri and quad faces.
1664 const int kNfaces = PyrGeom::kNfaces;
1665 const int kNtfaces = PyrGeom::kNtfaces;
1666 const int kNqfaces = PyrGeom::kNqfaces;
1667 Geometry2DSharedPtr faces[kNfaces];
1668 int Nfaces = 0;
1669 int Ntfaces = 0;
1670 int Nqfaces = 0;
1671
1672 /// Fill the arrays and make sure there aren't too many
1673 /// faces.
1674 std::stringstream errorstring;
1675 errorstring << "Element " << indx << " must have " << kNtfaces
1676 << " triangle face(s), and " << kNqfaces
1677 << " quadrilateral face(s).";
1678 for (int i = 0; i < kNfaces; i++)
1679 {
1680 int faceID;
1681 elementDataStrm >> faceID;
1682 Geometry2DSharedPtr face =
1683 m_meshGraph->GetGeometry2D(faceID);
1684 if (face == Geometry2DSharedPtr() ||
1685 (face->GetShapeType() != LibUtilities::eTriangle &&
1686 face->GetShapeType() != LibUtilities::eQuadrilateral))
1687 {
1688 std::stringstream errorstring;
1689 errorstring << "Element " << indx
1690 << " has invalid face: " << faceID;
1691 ASSERTL0(false, errorstring.str().c_str());
1692 }
1693 else if (face->GetShapeType() == LibUtilities::eTriangle)
1694 {
1695 ASSERTL0(Ntfaces < kNtfaces, errorstring.str().c_str());
1696 faces[Nfaces++] = static_pointer_cast<TriGeom>(face);
1697 Ntfaces++;
1698 }
1699 else if (face->GetShapeType() ==
1701 {
1702 ASSERTL0(Nqfaces < kNqfaces, errorstring.str().c_str());
1703 faces[Nfaces++] = static_pointer_cast<QuadGeom>(face);
1704 Nqfaces++;
1705 }
1706 }
1707
1708 /// Make sure all of the face indicies could be read, and
1709 /// that there weren't too few.
1710 ASSERTL0(
1711 !elementDataStrm.fail(),
1712 (std::string("Unable to read element data for PYRAMID: ") +
1713 elementStr)
1714 .c_str());
1715 ASSERTL0(Ntfaces == kNtfaces, errorstring.str().c_str());
1716 ASSERTL0(Nqfaces == kNqfaces, errorstring.str().c_str());
1717
1718 PyrGeomSharedPtr pyrgeom(
1720
1721 pyrGeoms[indx] = pyrgeom;
1722 m_meshGraph->PopulateFaceToElMap(pyrgeom, kNfaces);
1723 }
1724 catch (...)
1725 {
1726 NEKERROR(
1728 (std::string("Unable to read element data for PYRAMID: ") +
1729 elementStr)
1730 .c_str());
1731 }
1732 }
1733 // Prism
1734 else if (elementType == "R")
1735 {
1736 try
1737 {
1738 /// Create arrays for the tri and quad faces.
1739 const int kNfaces = PrismGeom::kNfaces;
1740 const int kNtfaces = PrismGeom::kNtfaces;
1741 const int kNqfaces = PrismGeom::kNqfaces;
1742 Geometry2DSharedPtr faces[kNfaces];
1743 int Ntfaces = 0;
1744 int Nqfaces = 0;
1745 int Nfaces = 0;
1746
1747 /// Fill the arrays and make sure there aren't too many
1748 /// faces.
1749 std::stringstream errorstring;
1750 errorstring << "Element " << indx << " must have " << kNtfaces
1751 << " triangle face(s), and " << kNqfaces
1752 << " quadrilateral face(s).";
1753
1754 for (int i = 0; i < kNfaces; i++)
1755 {
1756 int faceID;
1757 elementDataStrm >> faceID;
1758 Geometry2DSharedPtr face =
1759 m_meshGraph->GetGeometry2D(faceID);
1760 if (face == Geometry2DSharedPtr() ||
1761 (face->GetShapeType() != LibUtilities::eTriangle &&
1762 face->GetShapeType() != LibUtilities::eQuadrilateral))
1763 {
1764 std::stringstream errorstring;
1765 errorstring << "Element " << indx
1766 << " has invalid face: " << faceID;
1767 ASSERTL0(false, errorstring.str().c_str());
1768 }
1769 else if (face->GetShapeType() == LibUtilities::eTriangle)
1770 {
1771 ASSERTL0(Ntfaces < kNtfaces, errorstring.str().c_str());
1772 faces[Nfaces++] =
1773 std::static_pointer_cast<TriGeom>(face);
1774 Ntfaces++;
1775 }
1776 else if (face->GetShapeType() ==
1778 {
1779 ASSERTL0(Nqfaces < kNqfaces, errorstring.str().c_str());
1780 faces[Nfaces++] =
1781 std::static_pointer_cast<QuadGeom>(face);
1782 Nqfaces++;
1783 }
1784 }
1785
1786 /// Make sure all of the face indices could be read, and
1787 /// that there weren't too few.
1788 ASSERTL0(
1789 !elementDataStrm.fail(),
1790 (std::string("Unable to read element data for PRISM: ") +
1791 elementStr)
1792 .c_str());
1793 ASSERTL0(Ntfaces == kNtfaces, errorstring.str().c_str());
1794 ASSERTL0(Nqfaces == kNqfaces, errorstring.str().c_str());
1795
1796 PrismGeomSharedPtr prismgeom(
1798
1799 prismGeoms[indx] = prismgeom;
1800 m_meshGraph->PopulateFaceToElMap(prismgeom, kNfaces);
1801 }
1802 catch (...)
1803 {
1804 NEKERROR(
1806 (std::string("Unable to read element data for PRISM: ") +
1807 elementStr)
1808 .c_str());
1809 }
1810 }
1811 // Hexahedral
1812 else if (elementType == "H")
1813 {
1814 try
1815 {
1816 /// Create arrays for the tri and quad faces.
1817 const int kNfaces = HexGeom::kNfaces;
1818 const int kNtfaces = HexGeom::kNtfaces;
1819 const int kNqfaces = HexGeom::kNqfaces;
1820 // TriGeomSharedPtr tfaces[kNtfaces];
1821 QuadGeomSharedPtr qfaces[kNqfaces];
1822 int Ntfaces = 0;
1823 int Nqfaces = 0;
1824
1825 /// Fill the arrays and make sure there aren't too many
1826 /// faces.
1827 std::stringstream errorstring;
1828 errorstring << "Element " << indx << " must have " << kNtfaces
1829 << " triangle face(s), and " << kNqfaces
1830 << " quadrilateral face(s).";
1831 for (int i = 0; i < kNfaces; i++)
1832 {
1833 int faceID;
1834 elementDataStrm >> faceID;
1835 Geometry2DSharedPtr face =
1836 m_meshGraph->GetGeometry2D(faceID);
1837 if (face == Geometry2DSharedPtr() ||
1838 (face->GetShapeType() != LibUtilities::eTriangle &&
1839 face->GetShapeType() != LibUtilities::eQuadrilateral))
1840 {
1841 std::stringstream errorstring;
1842 errorstring << "Element " << indx
1843 << " has invalid face: " << faceID;
1844 ASSERTL0(false, errorstring.str().c_str());
1845 }
1846 else if (face->GetShapeType() == LibUtilities::eTriangle)
1847 {
1848 ASSERTL0(Ntfaces < kNtfaces, errorstring.str().c_str());
1849 // tfaces[Ntfaces++] =
1850 // boost::static_pointer_cast<TriGeom>(face);
1851 }
1852 else if (face->GetShapeType() ==
1854 {
1855 ASSERTL0(Nqfaces < kNqfaces, errorstring.str().c_str());
1856 qfaces[Nqfaces++] =
1857 std::static_pointer_cast<QuadGeom>(face);
1858 }
1859 }
1860
1861 /// Make sure all of the face indicies could be read, and
1862 /// that there weren't too few.
1863 ASSERTL0(!elementDataStrm.fail(),
1864 (std::string(
1865 "Unable to read element data for HEXAHEDRAL: ") +
1866 elementStr)
1867 .c_str());
1868 ASSERTL0(Ntfaces == kNtfaces, errorstring.str().c_str());
1869 ASSERTL0(Nqfaces == kNqfaces, errorstring.str().c_str());
1870
1871 HexGeomSharedPtr hexgeom(
1873
1874 hexGeoms[indx] = hexgeom;
1875 m_meshGraph->PopulateFaceToElMap(hexgeom, kNfaces);
1876 }
1877 catch (...)
1878 {
1880 (std::string(
1881 "Unable to read element data for HEXAHEDRAL: ") +
1882 elementStr)
1883 .c_str());
1884 }
1885 }
1886 /// Keep looking
1887 element = element->NextSiblingElement();
1888 }
1889}
1890
1892{
1893 auto &meshComposites = m_meshGraph->GetComposites();
1894 auto &compositesLabels = m_meshGraph->GetCompositesLabels();
1895
1896 TiXmlElement *field = nullptr;
1897
1898 /// Look for elements in ELEMENT block.
1899 field = m_xmlGeom->FirstChildElement("COMPOSITE");
1900
1901 ASSERTL0(field, "Unable to find COMPOSITE tag in file.");
1902
1903 TiXmlElement *node = field->FirstChildElement("C");
1904
1905 // Sequential counter for the composite numbers.
1906 int nextCompositeNumber = -1;
1907
1908 while (node)
1909 {
1910 /// All elements are of the form: "<? ID="#"> ... </?>", with
1911 /// ? being the element type.
1912
1913 nextCompositeNumber++;
1914
1915 int indx;
1916 int err = node->QueryIntAttribute("ID", &indx);
1917 ASSERTL0(err == TIXML_SUCCESS, "Unable to read attribute ID.");
1918 // ASSERTL0(indx == nextCompositeNumber, "Composite IDs must begin with
1919 // zero and be sequential.");
1920
1921 TiXmlNode *compositeChild = node->FirstChild();
1922 // This is primarily to skip comments that may be present.
1923 // Comments appear as nodes just like elements.
1924 // We are specifically looking for text in the body
1925 // of the definition.
1926 while (compositeChild &&
1927 compositeChild->Type() != TiXmlNode::TINYXML_TEXT)
1928 {
1929 compositeChild = compositeChild->NextSibling();
1930 }
1931
1932 ASSERTL0(compositeChild, "Unable to read composite definition body.");
1933 std::string compositeStr = compositeChild->ToText()->ValueStr();
1934
1935 /// Parse out the element components corresponding to type of element.
1936
1937 std::istringstream compositeDataStrm(compositeStr.c_str());
1938
1939 try
1940 {
1941 bool first = true;
1942 std::string prevCompositeElementStr;
1943
1944 while (!compositeDataStrm.fail())
1945 {
1946 std::string compositeElementStr;
1947 compositeDataStrm >> compositeElementStr;
1948
1949 if (!compositeDataStrm.fail())
1950 {
1951 if (first)
1952 {
1953 first = false;
1954 meshComposites[indx] =
1956 }
1957
1958 if (compositeElementStr.length() > 0)
1959 {
1960 ResolveGeomRef(prevCompositeElementStr,
1961 compositeElementStr,
1962 meshComposites[indx]);
1963 }
1964 prevCompositeElementStr = compositeElementStr;
1965 }
1966 }
1967 }
1968 catch (...)
1969 {
1970 NEKERROR(
1972 (std::string("Unable to read COMPOSITE data for composite: ") +
1973 compositeStr)
1974 .c_str());
1975 }
1976
1977 // Read optional name as string and save to compositeLabels if exists
1978 std::string name;
1979 err = node->QueryStringAttribute("NAME", &name);
1980 if (err == TIXML_SUCCESS)
1981 {
1982 compositesLabels[indx] = name;
1983 }
1984
1985 /// Keep looking for additional composite definitions.
1986 node = node->NextSiblingElement("C");
1987 }
1988
1989 ASSERTL0(nextCompositeNumber >= 0,
1990 "At least one composite must be specified.");
1991}
1992
1993void MeshGraphIOXml::ResolveGeomRef(const std::string &prevToken,
1994 const std::string &token,
1995 CompositeSharedPtr &composite)
1996{
1997 int meshDimension = m_meshGraph->GetMeshDimension();
1998
1999 switch (meshDimension)
2000 {
2001 case 1:
2002 ResolveGeomRef1D(prevToken, token, composite);
2003 break;
2004 case 2:
2005 ResolveGeomRef2D(prevToken, token, composite);
2006 break;
2007 case 3:
2008 ResolveGeomRef3D(prevToken, token, composite);
2009 break;
2010 }
2011}
2012
2013void MeshGraphIOXml::ResolveGeomRef1D(const std::string &prevToken,
2014 const std::string &token,
2015 CompositeSharedPtr &composite)
2016{
2017 auto &vertSet = m_meshGraph->GetAllPointGeoms();
2018 auto &segGeoms = m_meshGraph->GetAllSegGeoms();
2019
2020 try
2021 {
2022 std::istringstream tokenStream(token);
2023 std::istringstream prevTokenStream(prevToken);
2024
2025 char type;
2026 char prevType;
2027
2028 tokenStream >> type;
2029
2030 std::string::size_type indxBeg = token.find_first_of('[') + 1;
2031 std::string::size_type indxEnd = token.find_last_of(']') - 1;
2032
2033 ASSERTL0(
2034 indxBeg <= indxEnd,
2035 (std::string("Error reading index definition:") + token).c_str());
2036
2037 std::string indxStr = token.substr(indxBeg, indxEnd - indxBeg + 1);
2038
2039 typedef vector<unsigned int> SeqVectorType;
2040 SeqVectorType seqVector;
2041
2042 if (!ParseUtils::GenerateSeqVector(indxStr.c_str(), seqVector))
2043 {
2045 (std::string("Ill-formed sequence definition: ") + indxStr)
2046 .c_str());
2047 }
2048
2049 prevTokenStream >> prevType;
2050
2051 // All composites must be of the same dimension.
2052 bool validSequence =
2053 (prevToken.empty() || // No previous, then current is just fine.
2054 (type == 'V' && prevType == 'V') ||
2055 (type == 'S' && prevType == 'S'));
2056
2057 ASSERTL0(validSequence,
2058 std::string("Invalid combination of composite items: ") +
2059 type + " and " + prevType + ".");
2060
2061 switch (type)
2062 {
2063 case 'V': // Vertex
2064 for (SeqVectorType::iterator iter = seqVector.begin();
2065 iter != seqVector.end(); ++iter)
2066 {
2067 if (vertSet.find(*iter) == vertSet.end())
2068 {
2070 "Unknown vertex index: " +
2071 std::to_string(*iter));
2072 }
2073 else
2074 {
2075 composite->m_geomVec.push_back(vertSet[*iter]);
2076 }
2077 }
2078 break;
2079
2080 case 'S': // Segment
2081 for (SeqVectorType::iterator iter = seqVector.begin();
2082 iter != seqVector.end(); ++iter)
2083 {
2084 if (segGeoms.find(*iter) == segGeoms.end())
2085 {
2087 "Unknown segment index: " +
2088 std::to_string(*iter));
2089 }
2090 else
2091 {
2092 composite->m_geomVec.push_back(segGeoms[*iter]);
2093 }
2094 }
2095 break;
2096
2097 default:
2099 "Unrecognized composite token: " + token);
2100 }
2101 }
2102 catch (...)
2103 {
2105 "Problem processing composite token: " + token);
2106 }
2107
2108 return;
2109}
2110
2111void MeshGraphIOXml::ResolveGeomRef2D(const std::string &prevToken,
2112 const std::string &token,
2113 CompositeSharedPtr &composite)
2114{
2115 PointGeomMap &vertSet = m_meshGraph->GetAllPointGeoms();
2116 SegGeomMap &segGeoms = m_meshGraph->GetAllSegGeoms();
2117 TriGeomMap &triGeoms = m_meshGraph->GetAllTriGeoms();
2118 QuadGeomMap &quadGeoms = m_meshGraph->GetAllQuadGeoms();
2119
2120 try
2121 {
2122 std::istringstream tokenStream(token);
2123 std::istringstream prevTokenStream(prevToken);
2124
2125 char type;
2126 char prevType;
2127
2128 tokenStream >> type;
2129
2130 std::string::size_type indxBeg = token.find_first_of('[') + 1;
2131 std::string::size_type indxEnd = token.find_last_of(']') - 1;
2132
2133 ASSERTL0(
2134 indxBeg <= indxEnd,
2135 (std::string("Error reading index definition:") + token).c_str());
2136
2137 std::string indxStr = token.substr(indxBeg, indxEnd - indxBeg + 1);
2138 std::vector<unsigned int> seqVector;
2139 std::vector<unsigned int>::iterator seqIter;
2140
2141 bool err = ParseUtils::GenerateSeqVector(indxStr.c_str(), seqVector);
2142
2143 ASSERTL0(err,
2144 (std::string("Error reading composite elements: ") + indxStr)
2145 .c_str());
2146
2147 prevTokenStream >> prevType;
2148
2149 // All composites must be of the same dimension.
2150 bool validSequence =
2151 (prevToken.empty() || // No previous, then current is just fine.
2152 (type == 'V' && prevType == 'V') ||
2153 (type == 'E' && prevType == 'E') ||
2154 ((type == 'T' || type == 'Q') &&
2155 (prevType == 'T' || prevType == 'Q')));
2156
2157 ASSERTL0(validSequence,
2158 std::string("Invalid combination of composite items: ") +
2159 type + " and " + prevType + ".");
2160
2161 switch (type)
2162 {
2163 case 'E': // Edge
2164 for (seqIter = seqVector.begin(); seqIter != seqVector.end();
2165 ++seqIter)
2166 {
2167 if (segGeoms.find(*seqIter) == segGeoms.end())
2168 {
2170 "Unknown edge index: " +
2171 std::to_string(*seqIter));
2172 }
2173 else
2174 {
2175 composite->m_geomVec.push_back(segGeoms[*seqIter]);
2176 }
2177 }
2178 break;
2179
2180 case 'T': // Triangle
2181 {
2182 for (seqIter = seqVector.begin(); seqIter != seqVector.end();
2183 ++seqIter)
2184 {
2185 if (triGeoms.count(*seqIter) == 0)
2186 {
2188 "Unknown triangle index: " +
2189 std::to_string(*seqIter));
2190 }
2191 else
2192 {
2193 if (m_meshGraph->CheckRange(*triGeoms[*seqIter]))
2194 {
2195 composite->m_geomVec.push_back(triGeoms[*seqIter]);
2196 }
2197 }
2198 }
2199 }
2200 break;
2201
2202 case 'Q': // Quad
2203 {
2204 for (seqIter = seqVector.begin(); seqIter != seqVector.end();
2205 ++seqIter)
2206 {
2207 if (quadGeoms.count(*seqIter) == 0)
2208 {
2209 NEKERROR(
2211 "Unknown quad index: " + std::to_string(*seqIter) +
2212 " in Composite section");
2213 }
2214 else
2215 {
2216 if (m_meshGraph->CheckRange(*quadGeoms[*seqIter]))
2217 {
2218 composite->m_geomVec.push_back(quadGeoms[*seqIter]);
2219 }
2220 }
2221 }
2222 }
2223 break;
2224
2225 case 'V': // Vertex
2226 for (seqIter = seqVector.begin(); seqIter != seqVector.end();
2227 ++seqIter)
2228 {
2229 if (*seqIter >= vertSet.size())
2230 {
2232 "Unknown vertex index: " +
2233 std::to_string(*seqIter));
2234 }
2235 else
2236 {
2237 composite->m_geomVec.push_back(vertSet[*seqIter]);
2238 }
2239 }
2240 break;
2241
2242 default:
2244 "Unrecognized composite token: " + token);
2245 }
2246 }
2247 catch (...)
2248 {
2250 "Problem processing composite token: " + token);
2251 }
2252
2253 return;
2254}
2255
2256void MeshGraphIOXml::ResolveGeomRef3D(const std::string &prevToken,
2257 const std::string &token,
2258 CompositeSharedPtr &composite)
2259{
2260 PointGeomMap &vertSet = m_meshGraph->GetAllPointGeoms();
2261 SegGeomMap &segGeoms = m_meshGraph->GetAllSegGeoms();
2262 TriGeomMap &triGeoms = m_meshGraph->GetAllTriGeoms();
2263 QuadGeomMap &quadGeoms = m_meshGraph->GetAllQuadGeoms();
2264 TetGeomMap &tetGeoms = m_meshGraph->GetAllTetGeoms();
2265 PyrGeomMap &pyrGeoms = m_meshGraph->GetAllPyrGeoms();
2266 PrismGeomMap &prismGeoms = m_meshGraph->GetAllPrismGeoms();
2267 HexGeomMap &hexGeoms = m_meshGraph->GetAllHexGeoms();
2268
2269 try
2270 {
2271 std::istringstream tokenStream(token);
2272 std::istringstream prevTokenStream(prevToken);
2273
2274 char type;
2275 char prevType;
2276
2277 tokenStream >> type;
2278
2279 std::string::size_type indxBeg = token.find_first_of('[') + 1;
2280 std::string::size_type indxEnd = token.find_last_of(']') - 1;
2281
2282 ASSERTL0(indxBeg <= indxEnd,
2283 "Error reading index definition: " + token);
2284
2285 std::string indxStr = token.substr(indxBeg, indxEnd - indxBeg + 1);
2286
2287 std::vector<unsigned int> seqVector;
2288 std::vector<unsigned int>::iterator seqIter;
2289
2290 bool err = ParseUtils::GenerateSeqVector(indxStr.c_str(), seqVector);
2291
2292 ASSERTL0(err, "Error reading composite elements: " + indxStr);
2293
2294 prevTokenStream >> prevType;
2295
2296 // All composites must be of the same dimension. This map makes things
2297 // clean to compare.
2298 map<char, int> typeMap;
2299 typeMap['V'] = 1; // Vertex
2300 typeMap['E'] = 1; // Edge
2301 typeMap['T'] = 2; // Triangle
2302 typeMap['Q'] = 2; // Quad
2303 typeMap['A'] = 3; // Tet
2304 typeMap['P'] = 3; // Pyramid
2305 typeMap['R'] = 3; // Prism
2306 typeMap['H'] = 3; // Hex
2307
2308 // Make sure only geoms of the same dimension are combined.
2309 bool validSequence =
2310 (prevToken.empty() || (typeMap[type] == typeMap[prevType]));
2311
2312 ASSERTL0(validSequence,
2313 std::string("Invalid combination of composite items: ") +
2314 type + " and " + prevType + ".");
2315
2316 switch (type)
2317 {
2318 case 'V': // Vertex
2319 for (seqIter = seqVector.begin(); seqIter != seqVector.end();
2320 ++seqIter)
2321 {
2322 if (vertSet.find(*seqIter) == vertSet.end())
2323 {
2325 "Unknown vertex index: " +
2326 std::to_string(*seqIter));
2327 }
2328 else
2329 {
2330 composite->m_geomVec.push_back(vertSet[*seqIter]);
2331 }
2332 }
2333 break;
2334
2335 case 'E': // Edge
2336 for (seqIter = seqVector.begin(); seqIter != seqVector.end();
2337 ++seqIter)
2338 {
2339 if (segGeoms.find(*seqIter) == segGeoms.end())
2340 {
2342 "Unknown edge index: " +
2343 std::to_string(*seqIter));
2344 }
2345 else
2346 {
2347 composite->m_geomVec.push_back(segGeoms[*seqIter]);
2348 }
2349 }
2350 break;
2351
2352 case 'F': // Face
2353 for (seqIter = seqVector.begin(); seqIter != seqVector.end();
2354 ++seqIter)
2355 {
2356 Geometry2DSharedPtr face =
2357 m_meshGraph->GetGeometry2D(*seqIter);
2358 if (face == Geometry2DSharedPtr())
2359 {
2361 "Unknown face index: " +
2362 std::to_string(*seqIter));
2363 }
2364 else
2365 {
2366 if (m_meshGraph->CheckRange(*face))
2367 {
2368 composite->m_geomVec.push_back(face);
2369 }
2370 }
2371 }
2372 break;
2373
2374 case 'T': // Triangle
2375 for (seqIter = seqVector.begin(); seqIter != seqVector.end();
2376 ++seqIter)
2377 {
2378 if (triGeoms.find(*seqIter) == triGeoms.end())
2379 {
2381 "Unknown triangle index: " +
2382 std::to_string(*seqIter));
2383 }
2384 else
2385 {
2386 if (m_meshGraph->CheckRange(*triGeoms[*seqIter]))
2387 {
2388 composite->m_geomVec.push_back(triGeoms[*seqIter]);
2389 }
2390 }
2391 }
2392 break;
2393
2394 case 'Q': // Quad
2395 for (seqIter = seqVector.begin(); seqIter != seqVector.end();
2396 ++seqIter)
2397 {
2398 if (quadGeoms.find(*seqIter) == quadGeoms.end())
2399 {
2401 "Unknown quad index: " +
2402 std::to_string(*seqIter));
2403 }
2404 else
2405 {
2406 if (m_meshGraph->CheckRange(*quadGeoms[*seqIter]))
2407 {
2408 composite->m_geomVec.push_back(quadGeoms[*seqIter]);
2409 }
2410 }
2411 }
2412 break;
2413
2414 // Tetrahedron
2415 case 'A':
2416 for (seqIter = seqVector.begin(); seqIter != seqVector.end();
2417 ++seqIter)
2418 {
2419 if (tetGeoms.find(*seqIter) == tetGeoms.end())
2420 {
2422 "Unknown tet index: " +
2423 std::to_string(*seqIter));
2424 }
2425 else
2426 {
2427 if (m_meshGraph->CheckRange(*tetGeoms[*seqIter]))
2428 {
2429 composite->m_geomVec.push_back(tetGeoms[*seqIter]);
2430 }
2431 }
2432 }
2433 break;
2434
2435 // Pyramid
2436 case 'P':
2437 for (seqIter = seqVector.begin(); seqIter != seqVector.end();
2438 ++seqIter)
2439 {
2440 if (pyrGeoms.find(*seqIter) == pyrGeoms.end())
2441 {
2443 "Unknown pyramid index: " +
2444 std::to_string(*seqIter));
2445 }
2446 else
2447 {
2448 if (m_meshGraph->CheckRange(*pyrGeoms[*seqIter]))
2449 {
2450 composite->m_geomVec.push_back(pyrGeoms[*seqIter]);
2451 }
2452 }
2453 }
2454 break;
2455
2456 // Prism
2457 case 'R':
2458 for (seqIter = seqVector.begin(); seqIter != seqVector.end();
2459 ++seqIter)
2460 {
2461 if (prismGeoms.find(*seqIter) == prismGeoms.end())
2462 {
2464 "Unknown prism index: " +
2465 std::to_string(*seqIter));
2466 }
2467 else
2468 {
2469 if (m_meshGraph->CheckRange(*prismGeoms[*seqIter]))
2470 {
2471 composite->m_geomVec.push_back(
2472 prismGeoms[*seqIter]);
2473 }
2474 }
2475 }
2476 break;
2477
2478 // Hex
2479 case 'H':
2480 for (seqIter = seqVector.begin(); seqIter != seqVector.end();
2481 ++seqIter)
2482 {
2483 if (hexGeoms.find(*seqIter) == hexGeoms.end())
2484 {
2486 "Unknown hex index: " +
2487 std::to_string(*seqIter));
2488 }
2489 else
2490 {
2491 if (m_meshGraph->CheckRange(*hexGeoms[*seqIter]))
2492 {
2493 composite->m_geomVec.push_back(hexGeoms[*seqIter]);
2494 }
2495 }
2496 }
2497 break;
2498
2499 default:
2501 "Unrecognized composite token: " + token);
2502 }
2503 }
2504 catch (...)
2505 {
2507 "Problem processing composite token: " + token);
2508 }
2509
2510 return;
2511}
2512
2513void MeshGraphIOXml::v_WriteVertices(TiXmlElement *geomTag, PointGeomMap &verts)
2514{
2515 TiXmlElement *vertTag = new TiXmlElement("VERTEX");
2516
2517 for (auto &i : verts)
2518 {
2519 stringstream s;
2520 s << scientific << setprecision(8) << (*i.second)(0) << " "
2521 << (*i.second)(1) << " " << (*i.second)(2);
2522 TiXmlElement *v = new TiXmlElement("V");
2523 v->SetAttribute("ID", i.second->GetGlobalID());
2524 v->LinkEndChild(new TiXmlText(s.str()));
2525 vertTag->LinkEndChild(v);
2526 }
2527
2528 geomTag->LinkEndChild(vertTag);
2529}
2530
2531void MeshGraphIOXml::v_WriteEdges(TiXmlElement *geomTag, SegGeomMap &edges)
2532{
2533 int meshDimension = m_meshGraph->GetMeshDimension();
2534
2535 TiXmlElement *edgeTag =
2536 new TiXmlElement(meshDimension == 1 ? "ELEMENT" : "EDGE");
2537 string tag = meshDimension == 1 ? "S" : "E";
2538
2539 for (auto &i : edges)
2540 {
2541 stringstream s;
2542 SegGeomSharedPtr seg = i.second;
2543 s << seg->GetVid(0) << " " << seg->GetVid(1);
2544 TiXmlElement *e = new TiXmlElement(tag);
2545 e->SetAttribute("ID", i.first);
2546 e->LinkEndChild(new TiXmlText(s.str()));
2547 edgeTag->LinkEndChild(e);
2548 }
2549
2550 geomTag->LinkEndChild(edgeTag);
2551}
2552
2553void MeshGraphIOXml::v_WriteTris(TiXmlElement *faceTag, TriGeomMap &tris)
2554{
2555 string tag = "T";
2556
2557 for (auto &i : tris)
2558 {
2559 stringstream s;
2560 TriGeomSharedPtr tri = i.second;
2561 s << tri->GetEid(0) << " " << tri->GetEid(1) << " " << tri->GetEid(2);
2562 TiXmlElement *t = new TiXmlElement(tag);
2563 t->SetAttribute("ID", i.first);
2564 t->LinkEndChild(new TiXmlText(s.str()));
2565 faceTag->LinkEndChild(t);
2566 }
2567}
2568
2569void MeshGraphIOXml::v_WriteQuads(TiXmlElement *faceTag, QuadGeomMap &quads)
2570{
2571 string tag = "Q";
2572
2573 for (auto &i : quads)
2574 {
2575 stringstream s;
2576 QuadGeomSharedPtr quad = i.second;
2577 s << quad->GetEid(0) << " " << quad->GetEid(1) << " " << quad->GetEid(2)
2578 << " " << quad->GetEid(3);
2579 TiXmlElement *q = new TiXmlElement(tag);
2580 q->SetAttribute("ID", i.first);
2581 q->LinkEndChild(new TiXmlText(s.str()));
2582 faceTag->LinkEndChild(q);
2583 }
2584}
2585
2586void MeshGraphIOXml::v_WriteHexs(TiXmlElement *elmtTag, HexGeomMap &hexs)
2587{
2588 string tag = "H";
2589
2590 for (auto &i : hexs)
2591 {
2592 stringstream s;
2593 HexGeomSharedPtr hex = i.second;
2594 s << hex->GetFid(0) << " " << hex->GetFid(1) << " " << hex->GetFid(2)
2595 << " " << hex->GetFid(3) << " " << hex->GetFid(4) << " "
2596 << hex->GetFid(5) << " ";
2597 TiXmlElement *h = new TiXmlElement(tag);
2598 h->SetAttribute("ID", i.first);
2599 h->LinkEndChild(new TiXmlText(s.str()));
2600 elmtTag->LinkEndChild(h);
2601 }
2602}
2603
2604void MeshGraphIOXml::v_WritePrisms(TiXmlElement *elmtTag, PrismGeomMap &pris)
2605{
2606 string tag = "R";
2607
2608 for (auto &i : pris)
2609 {
2610 stringstream s;
2611 PrismGeomSharedPtr prism = i.second;
2612 s << prism->GetFid(0) << " " << prism->GetFid(1) << " "
2613 << prism->GetFid(2) << " " << prism->GetFid(3) << " "
2614 << prism->GetFid(4) << " ";
2615 TiXmlElement *p = new TiXmlElement(tag);
2616 p->SetAttribute("ID", i.first);
2617 p->LinkEndChild(new TiXmlText(s.str()));
2618 elmtTag->LinkEndChild(p);
2619 }
2620}
2621
2622void MeshGraphIOXml::v_WritePyrs(TiXmlElement *elmtTag, PyrGeomMap &pyrs)
2623{
2624 string tag = "P";
2625
2626 for (auto &i : pyrs)
2627 {
2628 stringstream s;
2629 PyrGeomSharedPtr pyr = i.second;
2630 s << pyr->GetFid(0) << " " << pyr->GetFid(1) << " " << pyr->GetFid(2)
2631 << " " << pyr->GetFid(3) << " " << pyr->GetFid(4) << " ";
2632 TiXmlElement *p = new TiXmlElement(tag);
2633 p->SetAttribute("ID", i.first);
2634 p->LinkEndChild(new TiXmlText(s.str()));
2635 elmtTag->LinkEndChild(p);
2636 }
2637}
2638
2639void MeshGraphIOXml::v_WriteTets(TiXmlElement *elmtTag, TetGeomMap &tets)
2640{
2641 string tag = "A";
2642
2643 for (auto &i : tets)
2644 {
2645 stringstream s;
2646 TetGeomSharedPtr tet = i.second;
2647 s << tet->GetFid(0) << " " << tet->GetFid(1) << " " << tet->GetFid(2)
2648 << " " << tet->GetFid(3) << " ";
2649 TiXmlElement *t = new TiXmlElement(tag);
2650 t->SetAttribute("ID", i.first);
2651 t->LinkEndChild(new TiXmlText(s.str()));
2652 elmtTag->LinkEndChild(t);
2653 }
2654}
2655
2656void MeshGraphIOXml::v_WriteCurves(TiXmlElement *geomTag, CurveMap &edges,
2657 CurveMap &faces)
2658{
2659 TiXmlElement *curveTag = new TiXmlElement("CURVED");
2660 CurveMap::iterator curveIt;
2661 int curveId = 0;
2662
2663 for (curveIt = edges.begin(); curveIt != edges.end(); ++curveIt)
2664 {
2665 CurveSharedPtr curve = curveIt->second;
2666 TiXmlElement *c = new TiXmlElement("E");
2667 stringstream s;
2668 s.precision(8);
2669
2670 for (int j = 0; j < curve->m_points.size(); ++j)
2671 {
2672 SpatialDomains::PointGeomSharedPtr p = curve->m_points[j];
2673 s << scientific << (*p)(0) << " " << (*p)(1) << " " << (*p)(2)
2674 << " ";
2675 }
2676
2677 c->SetAttribute("ID", curveId++);
2678 c->SetAttribute("EDGEID", curve->m_curveID);
2679 c->SetAttribute("NUMPOINTS", curve->m_points.size());
2680 c->SetAttribute("TYPE", LibUtilities::kPointsTypeStr[curve->m_ptype]);
2681 c->LinkEndChild(new TiXmlText(s.str()));
2682 curveTag->LinkEndChild(c);
2683 }
2684
2685 for (curveIt = faces.begin(); curveIt != faces.end(); ++curveIt)
2686 {
2687 CurveSharedPtr curve = curveIt->second;
2688 TiXmlElement *c = new TiXmlElement("F");
2689 stringstream s;
2690 s.precision(8);
2691
2692 for (int j = 0; j < curve->m_points.size(); ++j)
2693 {
2694 SpatialDomains::PointGeomSharedPtr p = curve->m_points[j];
2695 s << scientific << (*p)(0) << " " << (*p)(1) << " " << (*p)(2)
2696 << " ";
2697 }
2698
2699 c->SetAttribute("ID", curveId++);
2700 c->SetAttribute("FACEID", curve->m_curveID);
2701 c->SetAttribute("NUMPOINTS", curve->m_points.size());
2702 c->SetAttribute("TYPE", LibUtilities::kPointsTypeStr[curve->m_ptype]);
2703 c->LinkEndChild(new TiXmlText(s.str()));
2704 curveTag->LinkEndChild(c);
2705 }
2706
2707 geomTag->LinkEndChild(curveTag);
2708}
2709
2710void MeshGraphIOXml::WriteComposites(TiXmlElement *geomTag, CompositeMap &comps,
2711 std::map<int, std::string> &compLabels)
2712{
2713 auto &compositesLabels = m_meshGraph->GetCompositesLabels();
2714 TiXmlElement *compTag = new TiXmlElement("COMPOSITE");
2715
2716 for (auto &cIt : comps)
2717 {
2718 if (cIt.second->m_geomVec.size() == 0)
2719 {
2720 continue;
2721 }
2722
2723 TiXmlElement *c = new TiXmlElement("C");
2724 c->SetAttribute("ID", cIt.first);
2725 if (!compositesLabels[cIt.first].empty())
2726 {
2727 c->SetAttribute("NAME", compLabels[cIt.first]);
2728 }
2729 c->LinkEndChild(new TiXmlText(GetCompositeString(cIt.second)));
2730 compTag->LinkEndChild(c);
2731 }
2732
2733 geomTag->LinkEndChild(compTag);
2734}
2735
2736void MeshGraphIOXml::WriteDomain(TiXmlElement *geomTag,
2737 std::map<int, CompositeMap> &domainMap)
2738{
2739 TiXmlElement *domTag = new TiXmlElement("DOMAIN");
2740
2741 vector<unsigned int> idxList;
2742 for (auto &domain : domainMap)
2743 {
2744 TiXmlElement *c = new TiXmlElement("D");
2745 idxList.clear();
2746 stringstream s;
2747 s << " "
2748 << "C"
2749 << "[";
2750
2751 for (const auto &elem : domain.second)
2752 {
2753 idxList.push_back(elem.first);
2754 }
2755
2756 s << ParseUtils::GenerateSeqString(idxList) << "] ";
2757 c->SetAttribute("ID", domain.first);
2758 c->LinkEndChild(new TiXmlText(s.str()));
2759 domTag->LinkEndChild(c);
2760 }
2761
2762 geomTag->LinkEndChild(domTag);
2763}
2764
2766{
2767 int meshDimension = m_meshGraph->GetMeshDimension();
2768 auto &meshComposites = m_meshGraph->GetComposites();
2769
2770 TiXmlElement *expTag = new TiXmlElement("EXPANSIONS");
2771
2772 for (auto it = meshComposites.begin(); it != meshComposites.end(); it++)
2773 {
2774 if (it->second->m_geomVec[0]->GetShapeDim() == meshDimension)
2775 {
2776 TiXmlElement *exp = new TiXmlElement("E");
2777 exp->SetAttribute("COMPOSITE",
2778 "C[" + boost::lexical_cast<string>(it->first) +
2779 "]");
2780 exp->SetAttribute("NUMMODES", 4);
2781 exp->SetAttribute("TYPE", "MODIFIED");
2782 exp->SetAttribute("FIELDS", "u");
2783
2784 expTag->LinkEndChild(exp);
2785 }
2786 }
2787 root->LinkEndChild(expTag);
2788}
2789
2790/**
2791 * @brief Write out an XML file containing the GEOMETRY block
2792 * representing this MeshGraph instance inside a NEKTAR tag.
2793 */
2795 const std::string &outfilename, bool defaultExp,
2796 const LibUtilities::FieldMetaDataMap &metadata)
2797{
2798 int meshDimension = m_meshGraph->GetMeshDimension();
2799 int spaceDimension = m_meshGraph->GetSpaceDimension();
2800
2801 // Create empty TinyXML document.
2802 TiXmlDocument doc;
2803 TiXmlDeclaration *decl = new TiXmlDeclaration("1.0", "utf-8", "");
2804 doc.LinkEndChild(decl);
2805
2806 TiXmlElement *root = new TiXmlElement("NEKTAR");
2807 doc.LinkEndChild(root);
2808 TiXmlElement *geomTag = new TiXmlElement("GEOMETRY");
2809 root->LinkEndChild(geomTag);
2810
2811 // Add provenance information using FieldIO library.
2813 new LibUtilities::XmlTagWriter(root)),
2814 metadata);
2815
2816 // Update attributes with dimensions.
2817 geomTag->SetAttribute("DIM", meshDimension);
2818 geomTag->SetAttribute("SPACE", spaceDimension);
2819
2820 if (m_session != nullptr && !m_session->GetComm()->IsSerial())
2821 {
2822 geomTag->SetAttribute("PARTITION", m_session->GetComm()->GetRank());
2823 }
2824
2825 // Clear existing elements.
2826 geomTag->Clear();
2827
2828 // Write out informatio
2829 v_WriteVertices(geomTag, m_meshGraph->GetAllPointGeoms());
2830 v_WriteEdges(geomTag, m_meshGraph->GetAllSegGeoms());
2831 if (meshDimension > 1)
2832 {
2833 TiXmlElement *faceTag =
2834 new TiXmlElement(meshDimension == 2 ? "ELEMENT" : "FACE");
2835
2836 v_WriteTris(faceTag, m_meshGraph->GetAllTriGeoms());
2837 v_WriteQuads(faceTag, m_meshGraph->GetAllQuadGeoms());
2838 geomTag->LinkEndChild(faceTag);
2839 }
2840 if (meshDimension > 2)
2841 {
2842 TiXmlElement *elmtTag = new TiXmlElement("ELEMENT");
2843
2844 v_WriteHexs(elmtTag, m_meshGraph->GetAllHexGeoms());
2845 v_WritePyrs(elmtTag, m_meshGraph->GetAllPyrGeoms());
2846 v_WritePrisms(elmtTag, m_meshGraph->GetAllPrismGeoms());
2847 v_WriteTets(elmtTag, m_meshGraph->GetAllTetGeoms());
2848
2849 geomTag->LinkEndChild(elmtTag);
2850 }
2851 v_WriteCurves(geomTag, m_meshGraph->GetCurvedEdges(),
2852 m_meshGraph->GetCurvedFaces());
2853 WriteComposites(geomTag, m_meshGraph->GetComposites(),
2854 m_meshGraph->GetCompositesLabels());
2855 WriteDomain(geomTag, m_meshGraph->GetDomain());
2856
2857 if (defaultExp)
2858 {
2860 }
2861
2862 auto &movement = m_meshGraph->GetMovement();
2863 if (movement)
2864 {
2865 movement->WriteMovement(root);
2866 }
2867
2868 // Save file.
2869 doc.SaveFile(outfilename);
2870}
2871
2872void MeshGraphIOXml::WriteXMLGeometry(std::string outname,
2873 vector<set<unsigned int>> elements,
2874 vector<unsigned int> partitions)
2875{
2876 // so in theory this function is used by the mesh partitioner
2877 // giving instructions on how to write out a paritioned mesh.
2878 // the theory goes that the elements stored in the meshgraph are the
2879 // "whole" mesh so based on the information from the elmements list
2880 // we can filter the mesh entities and write some individual files
2881 // hopefully
2882
2883 // this is xml so we are going to write a directory with lots of
2884 // xml files
2885 string dirname = outname + "_xml";
2886 fs::path pdirname(dirname);
2887
2888 if (!fs::is_directory(dirname))
2889 {
2890 fs::create_directory(dirname);
2891 }
2892
2893 ASSERTL0(elements.size() == partitions.size(),
2894 "error in partitioned information");
2895
2896 for (int i = 0; i < partitions.size(); i++)
2897 {
2898 TiXmlDocument doc;
2899 TiXmlDeclaration *decl = new TiXmlDeclaration("1.0", "utf-8", "");
2900 doc.LinkEndChild(decl);
2901
2902 TiXmlElement *root = doc.FirstChildElement("NEKTAR");
2903 TiXmlElement *geomTag;
2904
2905 // Try to find existing NEKTAR tag.
2906 if (!root)
2907 {
2908 root = new TiXmlElement("NEKTAR");
2909 doc.LinkEndChild(root);
2910
2911 geomTag = new TiXmlElement("GEOMETRY");
2912 root->LinkEndChild(geomTag);
2913 }
2914 else
2915 {
2916 // Try to find existing GEOMETRY tag.
2917 geomTag = root->FirstChildElement("GEOMETRY");
2918
2919 if (!geomTag)
2920 {
2921 geomTag = new TiXmlElement("GEOMETRY");
2922 root->LinkEndChild(geomTag);
2923 }
2924 }
2925
2926 int meshDimension = m_meshGraph->GetMeshDimension();
2927 int spaceDimension = m_meshGraph->GetSpaceDimension();
2928
2929 geomTag->SetAttribute("DIM", meshDimension);
2930 geomTag->SetAttribute("SPACE", spaceDimension);
2931 geomTag->SetAttribute("PARTITION", partitions[i]);
2932
2933 // Add Mesh //
2934 // Get the elements
2935 auto &vertSet = m_meshGraph->GetAllPointGeoms();
2936 auto &segGeoms = m_meshGraph->GetAllSegGeoms();
2937 auto &triGeoms = m_meshGraph->GetAllTriGeoms();
2938 auto &quadGeoms = m_meshGraph->GetAllQuadGeoms();
2939 auto &tetGeoms = m_meshGraph->GetAllTetGeoms();
2940 auto &pyrGeoms = m_meshGraph->GetAllPyrGeoms();
2941 auto &prismGeoms = m_meshGraph->GetAllPrismGeoms();
2942 auto &hexGeoms = m_meshGraph->GetAllHexGeoms();
2943 auto &curvedEdges = m_meshGraph->GetCurvedEdges();
2944 auto &curvedFaces = m_meshGraph->GetCurvedFaces();
2945 auto &meshComposites = m_meshGraph->GetComposites();
2946 auto &globalDomain = m_meshGraph->GetDomain();
2947
2948 HexGeomMap localHex;
2949 PyrGeomMap localPyr;
2950 PrismGeomMap localPrism;
2951 TetGeomMap localTet;
2952 TriGeomMap localTri;
2953 QuadGeomMap localQuad;
2954 SegGeomMap localEdge;
2955 PointGeomMap localVert;
2956 CurveMap localCurveEdge;
2957 CurveMap localCurveFace;
2958
2959 vector<set<unsigned int>> entityIds(4);
2960 entityIds[meshDimension] = elements[i];
2961
2962 switch (meshDimension)
2963 {
2964 case 3:
2965 {
2966 for (auto &j : entityIds[3])
2967 {
2969 if (hexGeoms.count(j))
2970 {
2971 g = hexGeoms[j];
2972 localHex[j] = hexGeoms[j];
2973 }
2974 else if (pyrGeoms.count(j))
2975 {
2976 g = pyrGeoms[j];
2977 localPyr[j] = pyrGeoms[j];
2978 }
2979 else if (prismGeoms.count(j))
2980 {
2981 g = prismGeoms[j];
2982 localPrism[j] = prismGeoms[j];
2983 }
2984 else if (tetGeoms.count(j))
2985 {
2986 g = tetGeoms[j];
2987 localTet[j] = tetGeoms[j];
2988 }
2989 else
2990 {
2991 ASSERTL0(false, "element in partition not found");
2992 }
2993
2994 for (int k = 0; k < g->GetNumFaces(); k++)
2995 {
2996 entityIds[2].insert(g->GetFid(k));
2997 }
2998 for (int k = 0; k < g->GetNumEdges(); k++)
2999 {
3000 entityIds[1].insert(g->GetEid(k));
3001 }
3002 for (int k = 0; k < g->GetNumVerts(); k++)
3003 {
3004 entityIds[0].insert(g->GetVid(k));
3005 }
3006 }
3007 }
3008 break;
3009 case 2:
3010 {
3011 for (auto &j : entityIds[2])
3012 {
3014 if (triGeoms.count(j))
3015 {
3016 g = triGeoms[j];
3017 localTri[j] = triGeoms[j];
3018 }
3019 else if (quadGeoms.count(j))
3020 {
3021 g = quadGeoms[j];
3022 localQuad[j] = quadGeoms[j];
3023 }
3024 else
3025 {
3026 ASSERTL0(false, "element in partition not found");
3027 }
3028
3029 for (int k = 0; k < g->GetNumEdges(); k++)
3030 {
3031 entityIds[1].insert(g->GetEid(k));
3032 }
3033 for (int k = 0; k < g->GetNumVerts(); k++)
3034 {
3035 entityIds[0].insert(g->GetVid(k));
3036 }
3037 }
3038 }
3039 break;
3040 case 1:
3041 {
3042 for (auto &j : entityIds[1])
3043 {
3045 if (segGeoms.count(j))
3046 {
3047 g = segGeoms[j];
3048 localEdge[j] = segGeoms[j];
3049 }
3050 else
3051 {
3052 ASSERTL0(false, "element in partition not found");
3053 }
3054
3055 for (int k = 0; k < g->GetNumVerts(); k++)
3056 {
3057 entityIds[0].insert(g->GetVid(k));
3058 }
3059 }
3060 }
3061 }
3062
3063 if (meshDimension > 2)
3064 {
3065 for (auto &j : entityIds[2])
3066 {
3067 if (triGeoms.count(j))
3068 {
3069 localTri[j] = triGeoms[j];
3070 }
3071 else if (quadGeoms.count(j))
3072 {
3073 localQuad[j] = quadGeoms[j];
3074 }
3075 else
3076 {
3077 ASSERTL0(false, "face not found");
3078 }
3079 }
3080 }
3081
3082 if (meshDimension > 1)
3083 {
3084 for (auto &j : entityIds[1])
3085 {
3086 if (segGeoms.count(j))
3087 {
3088 localEdge[j] = segGeoms[j];
3089 }
3090 else
3091 {
3092 ASSERTL0(false, "edge not found");
3093 }
3094 }
3095 }
3096
3097 for (auto &j : entityIds[0])
3098 {
3099 if (vertSet.count(j))
3100 {
3101 localVert[j] = vertSet[j];
3102 }
3103 else
3104 {
3105 ASSERTL0(false, "vert not found");
3106 }
3107 }
3108
3109 v_WriteVertices(geomTag, localVert);
3110 v_WriteEdges(geomTag, localEdge);
3111 if (meshDimension > 1)
3112 {
3113 TiXmlElement *faceTag =
3114 new TiXmlElement(meshDimension == 2 ? "ELEMENT" : "FACE");
3115
3116 v_WriteTris(faceTag, localTri);
3117 v_WriteQuads(faceTag, localQuad);
3118 geomTag->LinkEndChild(faceTag);
3119 }
3120 if (meshDimension > 2)
3121 {
3122 TiXmlElement *elmtTag = new TiXmlElement("ELEMENT");
3123
3124 v_WriteHexs(elmtTag, localHex);
3125 v_WritePyrs(elmtTag, localPyr);
3126 v_WritePrisms(elmtTag, localPrism);
3127 v_WriteTets(elmtTag, localTet);
3128
3129 geomTag->LinkEndChild(elmtTag);
3130 }
3131
3132 for (auto &j : localTri)
3133 {
3134 if (curvedFaces.count(j.first))
3135 {
3136 localCurveFace[j.first] = curvedFaces[j.first];
3137 }
3138 }
3139 for (auto &j : localQuad)
3140 {
3141 if (curvedFaces.count(j.first))
3142 {
3143 localCurveFace[j.first] = curvedFaces[j.first];
3144 }
3145 }
3146 for (auto &j : localEdge)
3147 {
3148 if (curvedEdges.count(j.first))
3149 {
3150 localCurveEdge[j.first] = curvedEdges[j.first];
3151 }
3152 }
3153
3154 v_WriteCurves(geomTag, localCurveEdge, localCurveFace);
3155
3156 CompositeMap localComp;
3157 std::map<int, std::string> localCompLabels;
3158
3159 for (auto &j : meshComposites)
3160 {
3162 int dim = j.second->m_geomVec[0]->GetShapeDim();
3163
3164 for (int k = 0; k < j.second->m_geomVec.size(); k++)
3165 {
3166 if (entityIds[dim].count(j.second->m_geomVec[k]->GetGlobalID()))
3167 {
3168 comp->m_geomVec.push_back(j.second->m_geomVec[k]);
3169 }
3170 }
3171
3172 if (comp->m_geomVec.size())
3173 {
3174 localComp[j.first] = comp;
3175 auto compositesLabels = m_meshGraph->GetCompositesLabels();
3176 if (!compositesLabels[j.first].empty())
3177 {
3178 localCompLabels[j.first] = compositesLabels[j.first];
3179 }
3180 }
3181 }
3182
3183 WriteComposites(geomTag, localComp, localCompLabels);
3184
3185 map<int, CompositeMap> domain;
3186 for (auto &j : localComp)
3187 {
3188 for (auto &dom : globalDomain)
3189 {
3190 for (auto &dIt : dom.second)
3191 {
3192 if (j.first == dIt.first)
3193 {
3194 domain[dom.first][j.first] = j.second;
3195 break;
3196 }
3197 }
3198 }
3199 }
3200
3201 WriteDomain(geomTag, domain);
3202
3203 if (m_session->DefinesElement("NEKTAR/CONDITIONS"))
3204 {
3205 std::set<int> vBndRegionIdList;
3206 TiXmlElement *vConditions =
3207 new TiXmlElement(*m_session->GetElement("Nektar/Conditions"));
3208 TiXmlElement *vBndRegions =
3209 vConditions->FirstChildElement("BOUNDARYREGIONS");
3210 // Use fine-level for mesh partition (Parallel-in-Time)
3212 TiXmlElement *vBndConditions =
3213 vConditions->FirstChildElement("BOUNDARYCONDITIONS");
3214 // Use fine-level for mesh partition (Parallel-in-Time)
3216 0);
3217 TiXmlElement *vItem;
3218
3219 if (vBndRegions)
3220 {
3221 // Check for parallel-in-time
3222 bool multiLevel =
3223 vConditions->FirstChildElement("BOUNDARYREGIONS")
3224 ->FirstChildElement("TIMELEVEL") != nullptr;
3225
3226 TiXmlElement *vNewBndRegions =
3227 multiLevel ? new TiXmlElement("TIMELEVEL")
3228 : new TiXmlElement("BOUNDARYREGIONS");
3229 vItem = vBndRegions->FirstChildElement();
3230 auto graph_bndRegOrder = m_meshGraph->GetBndRegionOrdering();
3231 while (vItem)
3232 {
3233 std::string vSeqStr =
3234 vItem->FirstChild()->ToText()->Value();
3235 std::string::size_type indxBeg =
3236 vSeqStr.find_first_of('[') + 1;
3237 std::string::size_type indxEnd =
3238 vSeqStr.find_last_of(']') - 1;
3239 vSeqStr = vSeqStr.substr(indxBeg, indxEnd - indxBeg + 1);
3240 std::vector<unsigned int> vSeq;
3241 ParseUtils::GenerateSeqVector(vSeqStr.c_str(), vSeq);
3242
3243 vector<unsigned int> idxList;
3244
3245 for (unsigned int i = 0; i < vSeq.size(); ++i)
3246 {
3247 if (localComp.find(vSeq[i]) != localComp.end())
3248 {
3249 idxList.push_back(vSeq[i]);
3250 }
3251 }
3252 int p = atoi(vItem->Attribute("ID"));
3253
3254 std::string vListStr =
3256
3257 if (vListStr.length() == 0)
3258 {
3259 TiXmlElement *tmp = vItem;
3260 vItem = vItem->NextSiblingElement();
3261 vBndRegions->RemoveChild(tmp);
3262 }
3263 else
3264 {
3265 vListStr = "C[" + vListStr + "]";
3266 TiXmlText *vList = new TiXmlText(vListStr);
3267 TiXmlElement *vNewElement = new TiXmlElement("B");
3268 vNewElement->SetAttribute("ID", p);
3269 vNewElement->LinkEndChild(vList);
3270 vNewBndRegions->LinkEndChild(vNewElement);
3271 vBndRegionIdList.insert(p);
3272 vItem = vItem->NextSiblingElement();
3273 }
3274
3275 // store original bnd region order
3276 m_bndRegOrder[p] = vSeq;
3277 graph_bndRegOrder[p] = vSeq;
3278 }
3279 if (multiLevel)
3280 {
3281 // Use fine-level for mesh partition (Parallel-in-Time)
3282 size_t timeLevel = 0;
3283 while (vBndRegions)
3284 {
3285 vNewBndRegions->SetAttribute("VALUE", timeLevel);
3286 vConditions->FirstChildElement("BOUNDARYREGIONS")
3287 ->ReplaceChild(vBndRegions, *vNewBndRegions);
3288 vBndRegions = vBndRegions->NextSiblingElement();
3289 timeLevel++;
3290 }
3291 }
3292 else
3293 {
3294 vConditions->ReplaceChild(vBndRegions, *vNewBndRegions);
3295 }
3296 }
3297
3298 if (vBndConditions)
3299 {
3300 vItem = vBndConditions->FirstChildElement();
3301 while (vItem)
3302 {
3303 std::set<int>::iterator x;
3304 if ((x = vBndRegionIdList.find(atoi(vItem->Attribute(
3305 "REF")))) != vBndRegionIdList.end())
3306 {
3307 vItem->SetAttribute("REF", *x);
3308 vItem = vItem->NextSiblingElement();
3309 }
3310 else
3311 {
3312 TiXmlElement *tmp = vItem;
3313 vItem = vItem->NextSiblingElement();
3314 vBndConditions->RemoveChild(tmp);
3315 }
3316 }
3317 }
3318 root->LinkEndChild(vConditions);
3319 }
3320
3321 // Distribute other sections of the XML to each process as is.
3322 TiXmlElement *vSrc =
3323 m_session->GetElement("Nektar")->FirstChildElement();
3324 while (vSrc)
3325 {
3326 std::string vName = boost::to_upper_copy(vSrc->ValueStr());
3327 if (vName != "GEOMETRY" && vName != "CONDITIONS")
3328 {
3329 root->LinkEndChild(new TiXmlElement(*vSrc));
3330 }
3331 vSrc = vSrc->NextSiblingElement();
3332 }
3333
3334 // Save Mesh
3335
3336 boost::format pad("P%1$07d.xml");
3337 pad % partitions[i];
3338 fs::path pFilename(pad.str());
3339
3340 fs::path fullpath = pdirname / pFilename;
3341 doc.SaveFile(LibUtilities::PortablePath(fullpath));
3342 }
3343}
3344
3346{
3347 auto &meshComposites = m_meshGraph->GetComposites();
3348 auto &domain = m_meshGraph->GetDomain();
3349
3351
3352 for (auto &c : meshComposites)
3353 {
3354 bool fillComp = true;
3355 for (auto &d : domain[0])
3356 {
3357 if (c.second == d.second)
3358 {
3359 fillComp = false;
3360 }
3361 }
3362 if (fillComp)
3363 {
3364 std::vector<unsigned int> ids;
3365 for (auto &elmt : c.second->m_geomVec)
3366 {
3367 ids.push_back(elmt->GetGlobalID());
3368 }
3369 ret[c.first] = ids;
3370 }
3371 }
3372
3373 return ret;
3374}
3375
3376} // namespace Nektar::SpatialDomains
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
static void AddInfoTag(TagWriterSharedPtr root, const FieldMetaDataMap &fieldmetadatamap)
Add provenance information to the field metadata map.
Definition: FieldIO.cpp:342
Interpreter class for the evaluation of mathematical expressions.
Definition: Interpreter.h:76
int DefineFunction(const std::string &vlist, const std::string &expr)
Defines a function for the purposes of evaluation.
NekDouble Evaluate(const int id)
Evaluate a function which depends only on constants and/or parameters.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
bool ModuleExists(tKey idKey)
Checks if a particular module is available.
static void GetXMLElementTimeLevel(TiXmlElement *&element, const size_t timeLevel, const bool enableCheck=true)
Get XML elment time level (Parallel-in-Time)
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
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
static const int kNqfaces
Definition: HexGeom.h:56
static const int kNtfaces
Definition: HexGeom.h:57
LibUtilities::SessionReaderSharedPtr m_session
Definition: MeshGraphIO.h:90
void ReadGeometry(LibUtilities::DomainRangeShPtr rng, bool fillGraph)
Definition: MeshGraphIO.h:74
std::string GetCompositeString(CompositeSharedPtr comp)
Returns a string representation of a composite.
CompositeDescriptor CreateCompositeDescriptor()
std::map< int, MeshEntity > CreateMeshEntities()
Create mesh entities for this graph.
virtual void v_WriteHexs(TiXmlElement *elmtTag, HexGeomMap &hexs)
void ResolveGeomRef2D(const std::string &prevToken, const std::string &token, CompositeSharedPtr &composite)
void ResolveGeomRef3D(const std::string &prevToken, const std::string &token, CompositeSharedPtr &composite)
void v_WriteGeometry(const std::string &outfilename, bool defaultExp=false, const LibUtilities::FieldMetaDataMap &metadata=LibUtilities::NullFieldMetaDataMap) override
Write out an XML file containing the GEOMETRY block representing this MeshGraph instance inside a NEK...
virtual void v_WriteCurves(TiXmlElement *geomTag, CurveMap &edges, CurveMap &faces)
void ResolveGeomRef(const std::string &prevToken, const std::string &token, CompositeSharedPtr &composite)
virtual void v_WritePyrs(TiXmlElement *elmtTag, PyrGeomMap &pyrs)
void v_PartitionMesh(LibUtilities::SessionReaderSharedPtr session) override
void WriteDefaultExpansion(TiXmlElement *root)
void v_ReadGeometry(LibUtilities::DomainRangeShPtr rng, bool fillGraph) override
virtual void v_WriteTets(TiXmlElement *elmtTag, TetGeomMap &tets)
void ResolveGeomRef1D(const std::string &prevToken, const std::string &token, CompositeSharedPtr &composite)
virtual void v_WriteEdges(TiXmlElement *geomTag, SegGeomMap &edges)
virtual void v_WritePrisms(TiXmlElement *elmtTag, PrismGeomMap &pris)
void WriteComposites(TiXmlElement *geomTag, CompositeMap &comps, std::map< int, std::string > &compLabels)
void WriteXMLGeometry(std::string outname, std::vector< std::set< unsigned int > > elements, std::vector< unsigned int > partitions)
virtual void v_WriteVertices(TiXmlElement *geomTag, PointGeomMap &verts)
virtual void v_WriteQuads(TiXmlElement *faceTag, QuadGeomMap &quads)
virtual void v_WriteTris(TiXmlElement *faceTag, TriGeomMap &tris)
void WriteDomain(TiXmlElement *geomTag, std::map< int, CompositeMap > &domain)
static MeshGraphIOSharedPtr create()
static const int kNfaces
Definition: PyrGeom.h:55
static const int kNtfaces
Definition: PyrGeom.h:54
static const int kNqfaces
Definition: PyrGeom.h:53
static const int kNedges
Definition: QuadGeom.h:74
static const int kNfaces
Definition: TetGeom.h:56
static const int kNqfaces
Definition: TetGeom.h:54
static const int kNtfaces
Definition: TetGeom.h:55
static const int kNedges
Get the orientation of face1.
Definition: TriGeom.h:70
std::map< std::string, std::string > FieldMetaDataMap
Definition: FieldIO.h:50
const std::string kPointsTypeStr[]
Definition: Foundations.hpp:52
static std::string PortablePath(const fs::path &path)
create portable path on different platforms for std::filesystem path.
Definition: Filesystem.hpp:56
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< DomainRange > DomainRangeShPtr
Definition: DomainRange.h:64
std::shared_ptr< XmlTagWriter > XmlTagWriterSharedPtr
Definition: FieldIOXml.h:166
static DomainRangeShPtr NullDomainRangeShPtr
Definition: DomainRange.h:65
@ SIZE_PointsType
Length of enum list.
Definition: PointsType.h:95
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
std::map< int, TriGeomSharedPtr > TriGeomMap
Definition: TriGeom.h:57
std::map< int, std::vector< unsigned int > > CompositeOrdering
Definition: MeshGraph.h:107
std::shared_ptr< QuadGeom > QuadGeomSharedPtr
Definition: HexGeom.h:45
std::map< int, PyrGeomSharedPtr > PyrGeomMap
Definition: PyrGeom.h:76
std::map< int, QuadGeomSharedPtr > QuadGeomMap
Definition: QuadGeom.h:52
std::shared_ptr< PrismGeom > PrismGeomSharedPtr
Definition: PrismGeom.h:82
MeshPartitionFactory & GetMeshPartitionFactory()
std::shared_ptr< Composite > CompositeSharedPtr
Definition: MeshGraph.h:135
std::map< int, SegGeomSharedPtr > SegGeomMap
Definition: SegGeom.h:49
std::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:58
std::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:59
std::shared_ptr< HexGeom > HexGeomSharedPtr
Definition: HexGeom.h:84
std::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:59
std::map< int, TetGeomSharedPtr > TetGeomMap
Definition: TetGeom.h:84
std::shared_ptr< PyrGeom > PyrGeomSharedPtr
Definition: PyrGeom.h:75
MeshGraphIOFactory & GetMeshGraphIOFactory()
Definition: MeshGraphIO.cpp:47
std::shared_ptr< TetGeom > TetGeomSharedPtr
Definition: TetGeom.h:82
std::map< int, PrismGeomSharedPtr > PrismGeomMap
Definition: PrismGeom.h:83
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< Geometry > GeometrySharedPtr
Definition: Geometry.h:51
std::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:56
std::map< int, HexGeomSharedPtr > HexGeomMap
Definition: HexGeom.h:85
std::map< int, PointGeomSharedPtr > PointGeomMap
Definition: PointGeom.h:51
std::map< int, CompositeSharedPtr > CompositeMap
Definition: MeshGraph.h:136
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:475
std::vector< double > d(NPUPPER *NPUPPER)
std::vector< double > q(NPUPPER *NPUPPER)
double NekDouble
STL namespace.