Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes | List of all members
Nektar::Utilities::InputSem Class Reference

#include <InputSem.h>

Inheritance diagram for Nektar::Utilities::InputSem:
[legend]

Public Member Functions

 InputSem (NekMeshUtils::MeshSharedPtr m)
 Initialises the InputSem class. More...
 
virtual ~InputSem ()
 
virtual void Process ()
 Process a Semtex session file. More...
 
- Public Member Functions inherited from Nektar::NekMeshUtils::InputModule
NEKMESHUTILS_EXPORT InputModule (MeshSharedPtr p_m)
 
NEKMESHUTILS_EXPORT void OpenStream ()
 Open a file for input. More...
 
- Public Member Functions inherited from Nektar::NekMeshUtils::Module
NEKMESHUTILS_EXPORT Module (MeshSharedPtr p_m)
 
NEKMESHUTILS_EXPORT void RegisterConfig (std::string key, std::string value=std::string())
 Register a configuration option with a module. More...
 
NEKMESHUTILS_EXPORT void PrintConfig ()
 Print out all configuration options for a module. More...
 
NEKMESHUTILS_EXPORT void SetDefaults ()
 Sets default configuration options for those which have not been set. More...
 
NEKMESHUTILS_EXPORT MeshSharedPtr GetMesh ()
 
virtual NEKMESHUTILS_EXPORT void ProcessVertices ()
 Extract element vertices. More...
 
virtual NEKMESHUTILS_EXPORT void ProcessEdges (bool ReprocessEdges=true)
 Extract element edges. More...
 
virtual NEKMESHUTILS_EXPORT void ProcessFaces (bool ReprocessFaces=true)
 Extract element faces. More...
 
virtual NEKMESHUTILS_EXPORT void ProcessElements ()
 Generate element IDs. More...
 
virtual NEKMESHUTILS_EXPORT void ProcessComposites ()
 Generate composites. More...
 
virtual NEKMESHUTILS_EXPORT void ClearElementLinks ()
 

Static Public Member Functions

static NekMeshUtils::ModuleSharedPtr create (NekMeshUtils::MeshSharedPtr m)
 Creates an instance of this class. More...
 

Static Public Attributes

static NekMeshUtils::ModuleKey className
 ModuleKey for class. More...
 

Private Member Functions

void insertEdge (int elmt, int side, int tagId)
 

Private Attributes

std::stringstream m_fileStream
 
std::map< std::string, std::streampos > sectionMap
 Maps Semtex sections to positions inside the input file. More...
 

Additional Inherited Members

- Protected Member Functions inherited from Nektar::NekMeshUtils::InputModule
NEKMESHUTILS_EXPORT void PrintSummary ()
 Print summary of elements. More...
 
- Protected Member Functions inherited from Nektar::NekMeshUtils::Module
NEKMESHUTILS_EXPORT void ReorderPrisms (PerMap &perFaces)
 Reorder node IDs so that prisms and tetrahedra are aligned correctly. More...
 
NEKMESHUTILS_EXPORT void PrismLines (int prism, PerMap &perFaces, std::set< int > &prismsDone, std::vector< ElementSharedPtr > &line)
 
- Protected Attributes inherited from Nektar::NekMeshUtils::InputModule
io::filtering_istream m_mshFile
 Input stream. More...
 
std::ifstream m_mshFileStream
 Input stream. More...
 
- Protected Attributes inherited from Nektar::NekMeshUtils::Module
MeshSharedPtr m_mesh
 Mesh object. More...
 
std::map< std::string, ConfigOptionm_config
 List of configuration values. More...
 

Detailed Description

Converter for Semtex session files.

Definition at line 48 of file InputSem.h.

Constructor & Destructor Documentation

◆ InputSem()

Nektar::Utilities::InputSem::InputSem ( NekMeshUtils::MeshSharedPtr  m)

Initialises the InputSem class.

Definition at line 55 of file InputSem.cpp.

55  : InputModule(m)
56 {
57 }
NEKMESHUTILS_EXPORT InputModule(MeshSharedPtr p_m)

◆ ~InputSem()

Nektar::Utilities::InputSem::~InputSem ( )
virtual

Definition at line 59 of file InputSem.cpp.

60 {
61 }

Member Function Documentation

◆ create()

static NekMeshUtils::ModuleSharedPtr Nektar::Utilities::InputSem::create ( NekMeshUtils::MeshSharedPtr  m)
inlinestatic

Creates an instance of this class.

Definition at line 56 of file InputSem.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

57  {
59  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

◆ insertEdge()

void Nektar::Utilities::InputSem::insertEdge ( int  elmt,
int  side,
int  tagId 
)
private

Definition at line 663 of file InputSem.cpp.

References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eSegment, Nektar::NekMeshUtils::GetElementFactory(), and Nektar::NekMeshUtils::Module::m_mesh.

Referenced by Process().

664 {
665  EdgeSharedPtr edge = m_mesh->m_element[2][elmt]->GetEdge(side);
666  vector<NodeSharedPtr> edgeNodes = edge->m_edgeNodes;
667  edgeNodes.insert(edgeNodes.begin(), edge->m_n2);
668  edgeNodes.insert(edgeNodes.begin(), edge->m_n1);
669  int order = edgeNodes.size() - 1;
670 
671  vector<int> tags;
672  tags.push_back(tagId);
673 
675  order,
676  order > 1,
677  false,
678  true,
681  LibUtilities::eSegment, conf, edgeNodes, tags);
682  m_mesh->m_element[1].push_back(E);
683 }
Basic information about an element.
Definition: ElementConfig.h:49
std::shared_ptr< Edge > EdgeSharedPtr
Shared pointer to an edge.
Definition: Edge.h:136
ElementFactory & GetElementFactory()
Definition: Element.cpp:44
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
std::shared_ptr< Element > ElementSharedPtr
Definition: Edge.h:49
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51

◆ Process()

void Nektar::Utilities::InputSem::Process ( )
virtual

Process a Semtex session file.

Semtex files are defined by a tokenized markup format. We first populate sectionMap which stores the location of the various tags in the session file so that they can be jumped to, since no ordering is defined. The converter only requires the NODES and ELEMENTS sections to exist, but can also read CURVES and SURFACES. High-order curves rely on the meshfile session.msh to be created with the Semtex utility meshpr first.

Parameters
pFilenameFilename of Semtex session to read.

Implements Nektar::NekMeshUtils::Module.

Definition at line 76 of file InputSem.cpp.

References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::SpatialDomains::eDirichlet, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::NekMeshUtils::eHOPCondition, Nektar::SpatialDomains::eNeumann, Nektar::SpatialDomains::ePeriodic, Nektar::LibUtilities::eQuadrilateral, Nektar::NekMeshUtils::GetElementFactory(), insertEdge(), Nektar::NekMeshUtils::Module::m_config, m_fileStream, Nektar::NekMeshUtils::Module::m_mesh, Nektar::NekMeshUtils::InputModule::m_mshFile, class_topology::Node, Nektar::NekMeshUtils::InputModule::OpenStream(), CellMLToNektar.cellml_metadata::p, Nektar::NekMeshUtils::InputModule::PrintSummary(), Nektar::NekMeshUtils::Module::ProcessComposites(), Nektar::NekMeshUtils::Module::ProcessEdges(), Nektar::NekMeshUtils::Module::ProcessElements(), Nektar::NekMeshUtils::Module::ProcessFaces(), Nektar::NekMeshUtils::Module::ProcessVertices(), and sectionMap.

77 {
78  // Open the file stream.
79  OpenStream();
80 
81  if (m_mesh->m_verbose)
82  {
83  cout << "InputSem: Start reading file..." << endl;
84  }
85 
86  // Read through input file and populate the section map.
87  string fileContents, line, word;
88  stringstream ss, ssFile;
89  streampos linePos;
90 
91  sectionMap["NODES"] = -1;
92  sectionMap["ELEMENTS"] = -1;
93  sectionMap["CURVES"] = -1;
94  sectionMap["SURFACES"] = -1;
95  sectionMap["GROUPS"] = -1;
96  sectionMap["BCS"] = -1;
97  sectionMap["FIELDS"] = -1;
98 
99  // We need to read entire file into a string and wrap around a stringstream,
100  // since boost::io::filtered_stream does not support seeking with zlib.
101  m_fileStream << m_mshFile.rdbuf();
102 
103  while (!m_fileStream.eof())
104  {
105  linePos = m_fileStream.tellg();
106  getline(m_fileStream, line);
107  ss.clear();
108  ss.str(line);
109  ss >> word;
110 
111  // Iterate over all tokens and see if section exists on this
112  // line.
113  for (auto &it : sectionMap)
114  {
115  if (word == "<" + it.first || word == "<" + it.first + ">")
116  {
117  sectionMap[it.first] = linePos;
118  }
119  }
120  }
121 
122  // Clear eofbit and go back to the beginning of the file.
123  m_fileStream.clear();
124  m_fileStream.seekg(0);
125 
126  // Check that required sections exist in the file.
127  if (sectionMap["NODES"] == std::streampos(-1))
128  {
129  cerr << "Unable to locate NODES section in session file." << endl;
130  abort();
131  }
132 
133  if (sectionMap["ELEMENTS"] == std::streampos(-1))
134  {
135  cerr << "Unable to locate ELEMENTS section in session file." << endl;
136  abort();
137  }
138 
139  if (sectionMap["SURFACES"] != std::streampos(-1))
140  {
141  if (sectionMap["BCS"] == std::streampos(-1))
142  {
143  cerr << "SURFACES section defined but BCS section not "
144  << "found." << endl;
145  abort();
146  }
147 
148  if (sectionMap["GROUPS"] == std::streampos(-1))
149  {
150  cerr << "SURFACES section defined but GROUPS section not "
151  << "found." << endl;
152  abort();
153  }
154 
155  if (sectionMap["FIELDS"] == std::streampos(-1))
156  {
157  cerr << "SURFACES section defined but FIELDS section not "
158  << "found." << endl;
159  abort();
160  }
161  }
162 
163  m_mesh->m_expDim = 0;
164  string tag;
165  int start, end, nVertices, nEntities, nCurves, nSurf, nGroups, nBCs;
166  int id, i, j, k;
167  vector<double> hoXData, hoYData;
169  ifstream homeshFile;
170 
171  // Begin by reading in list of nodes which define the linear
172  // elements.
173  m_fileStream.seekg(sectionMap["NODES"]);
174  getline(m_fileStream, line);
175  ss.clear();
176  ss.str(line);
177  ss >> word;
178 
179  tag = ss.str();
180  start = tag.find_first_of('=');
181  end = tag.find_first_of('>');
182  nVertices = atoi(tag.substr(start + 1, end).c_str());
183 
184  i = id = 0;
185  while (i < nVertices)
186  {
187  getline(m_fileStream, line);
188  if (line.length() < 7)
189  continue;
190  ss.clear();
191  ss.str(line);
192  double x = 0, y = 0, z = 0;
193  ss >> id >> x >> y >> z;
194 
195  if ((y * y) > 0.000001 && m_mesh->m_spaceDim != 3)
196  {
197  m_mesh->m_spaceDim = 2;
198  }
199  if ((z * z) > 0.000001)
200  {
201  m_mesh->m_spaceDim = 3;
202  }
203  id -= 1; // counter starts at 0
204  m_mesh->m_node.push_back(
205  std::shared_ptr<Node>(new Node(id, x, y, z)));
206  ++i;
207  }
208 
209  // Now read in elements
210  m_fileStream.seekg(sectionMap["ELEMENTS"]);
211  getline(m_fileStream, line);
212  ss.clear();
213  ss.str(line);
214  ss >> word;
215 
216  tag = ss.str();
217  start = tag.find_first_of('=');
218  end = tag.find_first_of('>');
219  nEntities = atoi(tag.substr(start + 1, end).c_str());
220 
221  i = id = 0;
222  while (i < nEntities)
223  {
224  getline(m_fileStream, line);
225  if (line.length() < 18)
226  {
227  continue;
228  }
229 
230  // Create element tags
231  vector<int> tags;
232  tags.push_back(0); // composite
233 
234  // Read element node list
235  ss.clear();
236  ss.str(line);
237  ss >> id >> word;
238  vector<NodeSharedPtr> nodeList;
239  for (j = 0; j < 4; ++j)
240  {
241  int node = 0;
242  ss >> node;
243  nodeList.push_back(m_mesh->m_node[node - 1]);
244  }
245 
246  // Create element
247  ElmtConfig conf(elType, 1, false, false);
248  ElementSharedPtr E =
249  GetElementFactory().CreateInstance(elType, conf, nodeList, tags);
250 
251  // Determine mesh expansion dimension
252  if (E->GetDim() > m_mesh->m_expDim)
253  {
254  m_mesh->m_expDim = E->GetDim();
255  }
256  m_mesh->m_element[E->GetDim()].push_back(E);
257  ++i;
258  }
259 
260  // Finally, process curves.
261  if (sectionMap["CURVES"] != std::streampos(-1))
262  {
263  int np, nel, nodeId = m_mesh->m_node.size();
264 
265  m_fileStream.seekg(sectionMap["CURVES"]);
266  getline(m_fileStream, line);
267  ss.clear();
268  ss.str(line);
269  ss >> word;
270 
271  tag = ss.str();
272  start = tag.find_first_of('=');
273  end = tag.find_first_of('>');
274  nCurves = atoi(tag.substr(start + 1, end).c_str());
275 
276  // Some session files have empty curves sections; if nCurves
277  // is 0, no nead to load high order mesh file.
278  if (nCurves > 0)
279  {
280  string fname = m_config["infile"].as<string>();
281  int ext = fname.find_last_of('.');
282  string meshfile = fname.substr(0, ext) + ".msh";
283 
284  homeshFile.open(meshfile.c_str());
285  if (!homeshFile.is_open())
286  {
287  cerr << "Cannot open or find mesh file: " << meshfile << endl
288  << "Make sure to run meshpr on your session "
289  << "file first." << endl;
290  abort();
291  }
292 
293  // Make sure we have matching header.
294  getline(homeshFile, line);
295  ss.clear();
296  ss.str(line);
297  ss >> np >> nel >> nel >> nel;
298 
299  if (nel != m_mesh->m_element[m_mesh->m_expDim].size())
300  {
301  cerr << "Number of elements mismatch in mesh file." << endl;
302  abort();
303  }
304 
305  // Now read in all mesh data. This is horribly inefficient
306  // since not all elements are curved, but it is the
307  // easiest way of finding element data.
308  hoXData.resize(nel * np * np);
309  hoYData.resize(nel * np * np);
310 
311  for (j = 0; j < nel * np * np; ++j)
312  {
313  getline(homeshFile, line);
314  ss.clear();
315  ss.str(line);
316  ss >> hoXData[j] >> hoYData[j];
317  }
318 
319  homeshFile.close();
320  }
321 
322  i = id = 0;
323  while (i < nCurves)
324  {
325  getline(m_fileStream, line);
326  if (line.length() < 18)
327  {
328  continue;
329  }
330  int elmt = 0, side = 0;
331  ss.clear();
332  ss.str(line);
333  ss >> id >> elmt >> side >> word;
334  id--;
335  elmt--;
336 
337  vector<NodeSharedPtr> edgeNodes;
338 
339  if (word != "<SPLINE>" && word != "<ARC>")
340  {
341  cerr << "Unknown curve tag: " << word << endl;
342  abort();
343  }
344 
345  // See if we have already retrieved high-order data
346  // for this elements; prevents unnecessary computation
347  // for elements with multiple curves.
348  if (m_mesh->m_element[2][elmt]->GetConf().m_order > 1)
349  {
350  ++i;
351  continue;
352  }
353 
354  // Now set high order data for requested element.
355  for (int side = 0; side < 4; ++side)
356  {
357  int offset = elmt * np * np;
358  int stride = 0;
359 
360  switch (side)
361  {
362  case 0: // Bottom edge
363  offset += 0;
364  stride = 1;
365  break;
366  case 1: // Right edge
367  offset += np - 1;
368  stride = np;
369  break;
370  case 2: // Top edge
371  offset += np * np - 1;
372  stride = -1;
373  break;
374  case 3: // Left edge
375  offset += np * (np - 1);
376  stride = -np;
377  break;
378  default:
379  cerr << "Unknown side for curve id " << id << endl;
380  abort();
381  }
382 
383  for (j = 1; j < np - 1; ++j, ++nodeId)
384  {
385  double x = hoXData[offset + j * stride];
386  double y = hoYData[offset + j * stride];
387  NodeSharedPtr n =
388  std::shared_ptr<Node>(new Node(nodeId, x, y, 0.0));
389  edgeNodes.push_back(n);
390  }
391  }
392 
393  // Add internal points.
394  for (j = 1; j < np - 1; ++j)
395  {
396  int offset = j * np + 1;
397  for (k = 1; k < np - 1; ++k, ++nodeId)
398  {
399  double x = hoXData[offset + k];
400  double y = hoYData[offset + k];
401  NodeSharedPtr n =
402  std::shared_ptr<Node>(new Node(nodeId, x, y, 0.0));
403  edgeNodes.push_back(n);
404  }
405  }
406 
407  // Grab existing element from list and retrieve tags and
408  // vertices; insert these into existing edge nodes.
409  ElementSharedPtr e = m_mesh->m_element[2][elmt];
410  vector<NodeSharedPtr> elvert = e->GetVertexList();
411  vector<int> tags = e->GetTagList();
412  edgeNodes.insert(edgeNodes.begin(), elvert.begin(), elvert.end());
413 
414  // Create new element and replace with an incomplete
415  // quadrilateral of the correct order.
416  ElmtConfig conf(elType,
417  np - 1,
418  true,
419  false,
420  true,
422  m_mesh->m_element[2][elmt] = GetElementFactory().CreateInstance(
423  elType, conf, edgeNodes, tags);
424 
425  ++i;
426  }
427  }
428 
429  // Process field names
430  if (sectionMap["FIELDS"] != std::streampos(-1))
431  {
432  m_fileStream.seekg(sectionMap["FIELDS"]);
433  getline(m_fileStream, line);
434  getline(m_fileStream, line);
435  ss.clear();
436  ss.str(line);
437 
438  while (ss >> tag)
439  {
440  m_mesh->m_fields.push_back(tag);
441  }
442  }
443 
444  // Process surfaces if they exist. This is deliberately done after
445  // curves to ensure high-order points are preserved.
446  if (sectionMap["SURFACES"] != std::streampos(-1))
447  {
448  map<string, int> conditionMap;
449  int maxTag = -1;
450 
451  // First read in list of groups, which defines each condition tag.
452  m_fileStream.seekg(sectionMap["GROUPS"]);
453  getline(m_fileStream, line);
454  ss.clear();
455  ss.str(line);
456  ss >> word;
457 
458  tag = ss.str();
459  start = tag.find_first_of('=');
460  end = tag.find_first_of('>');
461  nGroups = atoi(tag.substr(start + 1, end).c_str());
462 
463  i = id = 0;
464  while (i < nGroups)
465  {
466  getline(m_fileStream, line);
467  ss.clear();
468  ss.str(line);
469  ss >> id >> tag;
470  conditionMap[tag] = i++;
471  }
472 
473  maxTag = i;
474 
475  // Now read in actual values for boundary conditions from BCS
476  // section.
477  m_fileStream.seekg(sectionMap["BCS"]);
478  getline(m_fileStream, line);
479  ss.clear();
480  ss.str(line);
481  ss >> word;
482 
483  tag = ss.str();
484  start = tag.find_first_of('=');
485  end = tag.find_first_of('>');
486  nBCs = atoi(tag.substr(start + 1, end).c_str());
487 
488  i = id = 0;
489  while (i < nBCs)
490  {
491  int nF;
492  string tmp;
494  getline(m_fileStream, line);
495  ss.clear();
496  ss.str(line);
497  ss >> id >> tag >> nF;
498 
499  p = ConditionSharedPtr(new Condition());
500  m_mesh->m_condition[conditionMap[tag]] = p;
501 
502  // Read boundary condition.
503  j = 0;
504  while (j < nF)
505  {
506  getline(m_fileStream, line);
507  ss.clear();
508  ss.str(line);
509  ss >> tmp;
510 
511  // First string should be condition type.
512  if (tmp == "<D>")
513  {
514  p->type.push_back(eDirichlet);
515  }
516  else if (tmp == "<N>")
517  {
518  p->type.push_back(eNeumann);
519  }
520  else if (tmp == "<H>")
521  {
522  p->type.push_back(eHOPCondition);
523  p->value.push_back("0");
524  p->field.push_back("p");
525  ++j;
526  continue;
527  }
528  else
529  {
530  cerr << "Unsupported boundary condition type " << tmp
531  << endl;
532  abort();
533  }
534 
535  // Second string should be field.
536  ss >> tmp;
537  p->field.push_back(tmp);
538 
539  // Third string should be equals sign.
540  ss >> tmp;
541  if (tmp != "=")
542  {
543  cerr << "Couldn't read boundary condition type " << tag
544  << endl;
545  abort();
546  }
547 
548  // Fourth string should be value. CAUTION: Assumes
549  // expression is defined without any spaces in it!
550  ss >> tmp;
551  p->value.push_back(tmp);
552 
553  ++j;
554  }
555 
556  // Finally set composite for condition. In this case, all
557  // composites will be lines so there is one set per
558  // composite.
559  p->m_composite.push_back(conditionMap[tag] + 1);
560 
561  ++i;
562  }
563 
564  // Finally read surface information.
565  m_fileStream.seekg(sectionMap["SURFACES"]);
566  getline(m_fileStream, line);
567  ss.clear();
568  ss.str(line);
569  ss >> word;
570 
571  tag = ss.str();
572  start = tag.find_first_of('=');
573  end = tag.find_first_of('>');
574  nSurf = atoi(tag.substr(start + 1, end).c_str());
575 
576  i = id = 0;
577  int elmt, side;
578  int periodicTagId = -1;
579 
580  set<pair<int, int> > visitedPeriodic;
581 
582  while (i < nSurf)
583  {
584  getline(m_fileStream, line);
585  ss.clear();
586  ss.str(line);
587  ss >> id >> elmt >> side >> word;
588  elmt--;
589  side--;
590 
591  if (word == "<P>")
592  {
593  // If this is the first periodic boundary condition
594  // encountered, then set up m_mesh->m_condition with two
595  // periodic conditions.
596  if (periodicTagId == -1)
597  {
598  periodicTagId = maxTag;
600  ConditionSharedPtr out =
602  for (j = 0; j < m_mesh->m_fields.size(); ++j)
603  {
604  in->type.push_back(ePeriodic);
605  out->type.push_back(ePeriodic);
606  in->field.push_back(m_mesh->m_fields[j]);
607  out->field.push_back(m_mesh->m_fields[j]);
608  in->value.push_back("[" + boost::lexical_cast<string>(
609  periodicTagId + 1) +
610  "]");
611  out->value.push_back(
612  "[" + boost::lexical_cast<string>(periodicTagId) +
613  "]");
614  }
615  in->m_composite.push_back(periodicTagId + 1);
616  out->m_composite.push_back(periodicTagId + 2);
617  m_mesh->m_condition[periodicTagId] = in;
618  m_mesh->m_condition[periodicTagId + 1] = out;
619  }
620 
621  int elmtB, sideB;
622 
623  ss >> elmtB >> sideB;
624  elmtB--;
625  sideB--;
626 
627  pair<int, int> c1(elmt, side);
628  pair<int, int> c2(elmtB, sideB);
629 
630  if (visitedPeriodic.count(c1) == 0 &&
631  visitedPeriodic.count(c2) == 0)
632  {
633  visitedPeriodic.insert(make_pair(elmtB, sideB));
634  visitedPeriodic.insert(make_pair(elmt, side));
635  insertEdge(elmt, side, periodicTagId + 1);
636  insertEdge(elmtB, sideB, periodicTagId + 2);
637  }
638  }
639  else if (word == "<B>")
640  {
641  ss >> tag;
642  insertEdge(elmt, side, conditionMap[tag] + 1);
643  }
644  else
645  {
646  cerr << "Unrecognised or unsupported tag " << word << endl;
647  abort();
648  }
649  ++i;
650  }
651  }
652 
653  PrintSummary();
654 
655  // Process rest of mesh.
656  ProcessVertices();
657  ProcessEdges();
658  ProcessFaces();
659  ProcessElements();
661 }
Defines a boundary condition.
Definition: Mesh.h:74
Basic information about an element.
Definition: ElementConfig.h:49
io::filtering_istream m_mshFile
Input stream.
std::stringstream m_fileStream
Definition: InputSem.h:64
ElementFactory & GetElementFactory()
Definition: Element.cpp:44
std::shared_ptr< Node > NodeSharedPtr
Definition: CADVert.h:49
NEKMESHUTILS_EXPORT void PrintSummary()
Print summary of elements.
NEKMESHUTILS_EXPORT void OpenStream()
Open a file for input.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
virtual NEKMESHUTILS_EXPORT void ProcessFaces(bool ReprocessFaces=true)
Extract element faces.
std::shared_ptr< Element > ElementSharedPtr
Definition: Edge.h:49
virtual NEKMESHUTILS_EXPORT void ProcessElements()
Generate element IDs.
std::map< std::string, std::streampos > sectionMap
Maps Semtex sections to positions inside the input file.
Definition: InputSem.h:68
std::map< std::string, ConfigOption > m_config
List of configuration values.
virtual NEKMESHUTILS_EXPORT void ProcessVertices()
Extract element vertices.
void insertEdge(int elmt, int side, int tagId)
Definition: InputSem.cpp:663
virtual NEKMESHUTILS_EXPORT void ProcessEdges(bool ReprocessEdges=true)
Extract element edges.
std::shared_ptr< Condition > ConditionSharedPtr
Definition: Mesh.h:85
virtual NEKMESHUTILS_EXPORT void ProcessComposites()
Generate composites.
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51

Member Data Documentation

◆ className

ModuleKey Nektar::Utilities::InputSem::className
static
Initial value:

ModuleKey for class.

Definition at line 61 of file InputSem.h.

◆ m_fileStream

std::stringstream Nektar::Utilities::InputSem::m_fileStream
private

Definition at line 64 of file InputSem.h.

Referenced by Process().

◆ sectionMap

std::map<std::string, std::streampos> Nektar::Utilities::InputSem::sectionMap
private

Maps Semtex sections to positions inside the input file.

Definition at line 68 of file InputSem.h.

Referenced by Process().