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