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