Nektar++
MeshGraph2D.cpp
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////////////
2 // File: MeshGraph2D.cpp
3 //
4 // For more information, please see: http://www.nektar.info/
5 //
6 // The MIT License
7 //
8 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
9 // Department of Aeronautics, Imperial College London (UK), and Scientific
10 // Computing and Imaging Institute, University of Utah (USA).
11 //
12 // License for the specific language governing rights and limitations under
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 
37 #include <SpatialDomains/SegGeom.h>
38 #include <SpatialDomains/TriGeom.h>
40 #include <tinyxml.h>
41 
42 namespace Nektar
43 {
44  namespace SpatialDomains
45  {
47  {
48  }
49 
51  {
52  }
53 
55  const DomainRangeShPtr &rng)
56  : MeshGraph(pSession,rng)
57  {
58  ReadGeometry(pSession->GetDocument());
59  ReadExpansions(pSession->GetDocument());
60  }
61 
62  void MeshGraph2D::ReadGeometry(const std::string &infilename)
63  {
64  TiXmlDocument doc(infilename);
65  bool loadOkay = doc.LoadFile();
66 
67  std::stringstream errstr;
68  errstr << "Unable to load file: " << infilename << "\n";
69  errstr << doc.ErrorDesc() << " (Line " << doc.ErrorRow()
70  << ", Column " << doc.ErrorCol() << ")";
71  ASSERTL0(loadOkay, errstr.str());
72 
73  ReadGeometry(doc);
74  }
75 
76  // \brief Read segments (and general MeshGraph) given TiXmlDocument.
77  void MeshGraph2D::ReadGeometry(TiXmlDocument &doc)
78  {
79  // Read mesh first
81  TiXmlHandle docHandle(&doc);
82 
83  TiXmlElement* mesh = NULL;
84 
85  /// Look for all geometry related data in GEOMETRY block.
86  mesh = docHandle.FirstChildElement("NEKTAR").FirstChildElement("GEOMETRY").Element();
87 
88  ASSERTL0(mesh, "Unable to find GEOMETRY tag in file.");
89 
90  ReadCurves(doc);
91  ReadEdges(doc);
92  ReadElements(doc);
93  ReadComposites(doc);
94  ReadDomain(doc);
95  }
96 
97  void MeshGraph2D::ReadEdges(TiXmlDocument &doc)
98  {
99  /// We know we have it since we made it this far.
100  TiXmlHandle docHandle(&doc);
101  TiXmlElement* mesh = docHandle.FirstChildElement("NEKTAR").FirstChildElement("GEOMETRY").Element();
102  TiXmlElement* field = NULL;
103 
104  /// Look for elements in ELEMENT block.
105  field = mesh->FirstChildElement("EDGE");
106 
107  ASSERTL0(field, "Unable to find EDGE tag in file.");
108 
109  /// All elements are of the form: "<E ID="#"> ... </E>", with
110  /// ? being the element type.
111  /// Read the ID field first.
112  TiXmlElement *edge = field->FirstChildElement("E");
113 
114  /// Since all edge data is one big text block, we need to
115  /// accumulate all TINYXML_TEXT data and then parse it. This
116  /// approach effectively skips all comments or other node
117  /// types since we only care about the edge list. We
118  /// cannot handle missing edge numbers as we could with
119  /// missing element numbers due to the text block format.
120  std::string edgeStr;
121  int indx;
123 
124  while(edge)
125  {
126  int err = edge->QueryIntAttribute("ID",&indx);
127  ASSERTL0(err == TIXML_SUCCESS, "Unable to read edge attribute ID.");
128 
129  TiXmlNode *child = edge->FirstChild();
130  edgeStr.clear();
131  if (child->Type() == TiXmlNode::TINYXML_TEXT)
132  {
133  edgeStr += child->ToText()->ValueStr();
134  }
135 
136  /// Now parse out the edges, three fields at a time.
137  int vertex1, vertex2;
138  std::istringstream edgeDataStrm(edgeStr.c_str());
139 
140  try
141  {
142  while (!edgeDataStrm.fail())
143  {
144  edgeDataStrm >> vertex1 >> vertex2;
145 
146  // Must check after the read because we may be
147  // at the end and not know it. If we are at
148  // the end we will add a duplicate of the last
149  // entry if we don't check here.
150  if (!edgeDataStrm.fail())
151  {
152  PointGeomSharedPtr vertices[2] = {GetVertex(vertex1), GetVertex(vertex2)};
153 
154  SegGeomSharedPtr edge;
155 
156  it = m_curvedEdges.find(indx);
157 
158  if(it == m_curvedEdges.end())
159  {
161  edge->SetGlobalID(indx); // Set global mesh id
162  }
163  else
164  {
165  edge = MemoryManager<SegGeom>::AllocateSharedPtr(indx, m_spaceDimension, vertices, it->second);
166  edge->SetGlobalID(indx); //Set global mesh id
167  }
168 
169  m_segGeoms[indx] = edge;
170  }
171  }
172  }
173  catch(...)
174  {
175  NEKERROR(ErrorUtil::efatal, (std::string("Unable to read edge data: ") + edgeStr).c_str());
176  }
177 
178  edge = edge->NextSiblingElement("E");
179  }
180 
181  }
182 
183  void MeshGraph2D::ReadElements(TiXmlDocument &doc)
184  {
185  /// We know we have it since we made it this far.
186  TiXmlHandle docHandle(&doc);
187  TiXmlElement* mesh = docHandle.FirstChildElement("NEKTAR").FirstChildElement("GEOMETRY").Element();
188  TiXmlElement* field = NULL;
189 
190  /// Look for elements in ELEMENT block.
191  field = mesh->FirstChildElement("ELEMENT");
192 
193  ASSERTL0(field, "Unable to find ELEMENT tag in file.");
194 
195  // Set up curve map for curved elements on an embedded manifold.
197 
198  /// All elements are of the form: "<? ID="#"> ... </?>", with
199  /// ? being the element type.
200 
201  TiXmlElement *element = field->FirstChildElement();
202 
203  while (element)
204  {
205  std::string elementType(element->ValueStr());
206 
207  ASSERTL0(elementType == "Q" || elementType == "T",
208  (std::string("Unknown 2D element type: ") + elementType).c_str());
209 
210  /// Read id attribute.
211  int indx;
212  int err = element->QueryIntAttribute("ID", &indx);
213  ASSERTL0(err == TIXML_SUCCESS, "Unable to read element attribute ID.");
214 
215  it = m_curvedFaces.find(indx);
216 
217  /// Read text element description.
218  TiXmlNode* elementChild = element->FirstChild();
219  std::string elementStr;
220  while(elementChild)
221  {
222  if (elementChild->Type() == TiXmlNode::TINYXML_TEXT)
223  {
224  elementStr += elementChild->ToText()->ValueStr();
225  }
226  elementChild = elementChild->NextSibling();
227  }
228 
229  ASSERTL0(!elementStr.empty(), "Unable to read element description body.");
230 
231  /// Parse out the element components corresponding to type of element.
232  if (elementType == "T")
233  {
234  // Read three edge numbers
235  int edge1, edge2, edge3;
236  std::istringstream elementDataStrm(elementStr.c_str());
237 
238  try
239  {
240  elementDataStrm >> edge1;
241  elementDataStrm >> edge2;
242  elementDataStrm >> edge3;
243 
244  ASSERTL0(!elementDataStrm.fail(), (std::string("Unable to read element data for TRIANGLE: ") + elementStr).c_str());
245 
246  /// Create a TriGeom to hold the new definition.
248  {
249  GetSegGeom(edge1),
250  GetSegGeom(edge2),
251  GetSegGeom(edge3)
252  };
253 
255  {
256  SegGeom::GetEdgeOrientation(*edges[0], *edges[1]),
257  SegGeom::GetEdgeOrientation(*edges[1], *edges[2]),
258  SegGeom::GetEdgeOrientation(*edges[2], *edges[0])
259  };
260 
261  TriGeomSharedPtr trigeom;
262  if (it == m_curvedFaces.end())
263  {
264  trigeom = MemoryManager<TriGeom>
266  edges,
267  edgeorient);
268  }
269  else
270  {
271  trigeom = MemoryManager<TriGeom>
273  edges,
274  edgeorient,
275  it->second);
276  }
277  trigeom->SetGlobalID(indx);
278 
279  m_triGeoms[indx] = trigeom;
280  }
281  catch(...)
282  {
283  NEKERROR(ErrorUtil::efatal, (std::string("Unable to read element data for TRIANGLE: ") + elementStr).c_str());
284  }
285  }
286  else if (elementType == "Q")
287  {
288  // Read four edge numbers
289  int edge1, edge2, edge3, edge4;
290  std::istringstream elementDataStrm(elementStr.c_str());
291 
292  try
293  {
294  elementDataStrm >> edge1;
295  elementDataStrm >> edge2;
296  elementDataStrm >> edge3;
297  elementDataStrm >> edge4;
298 
299  ASSERTL0(!elementDataStrm.fail(), (std::string("Unable to read element data for QUAD: ") + elementStr).c_str());
300 
301  /// Create a QuadGeom to hold the new definition.
303  {GetSegGeom(edge1),GetSegGeom(edge2),
304  GetSegGeom(edge3),GetSegGeom(edge4)};
305 
307  {
308  SegGeom::GetEdgeOrientation(*edges[0], *edges[1]),
309  SegGeom::GetEdgeOrientation(*edges[1], *edges[2]),
310  SegGeom::GetEdgeOrientation(*edges[2], *edges[3]),
311  SegGeom::GetEdgeOrientation(*edges[3], *edges[0])
312  };
313 
314  QuadGeomSharedPtr quadgeom;
315  if (it == m_curvedFaces.end())
316  {
317  quadgeom = MemoryManager<QuadGeom>
319  edges,
320  edgeorient);
321  }
322  else
323  {
324  quadgeom = MemoryManager<QuadGeom>
326  edges,
327  edgeorient,
328  it->second);
329  }
330  quadgeom->SetGlobalID(indx);
331 
332  m_quadGeoms[indx] = quadgeom;
333 
334  }
335  catch(...)
336  {
337  NEKERROR(ErrorUtil::efatal,(std::string("Unable to read element data for QUAD: ") + elementStr).c_str());
338  }
339  }
340 
341  /// Keep looking
342  element = element->NextSiblingElement();
343  }
344  }
345 
346  void MeshGraph2D::ReadComposites(TiXmlDocument &doc)
347  {
348  TiXmlHandle docHandle(&doc);
349 
350  /// We know we have it since we made it this far.
351  TiXmlElement* mesh = docHandle.FirstChildElement("NEKTAR").FirstChildElement("GEOMETRY").Element();
352  TiXmlElement* field = NULL;
353 
354  ASSERTL0(mesh, "Unable to find GEOMETRY tag in file.");
355 
356  /// Look for elements in ELEMENT block.
357  field = mesh->FirstChildElement("COMPOSITE");
358 
359  ASSERTL0(field, "Unable to find COMPOSITE tag in file.");
360 
361  int nextCompositeNumber = -1;
362 
363  /// All elements are of the form: "<C ID = "N"> ... </C>".
364  /// Read the ID field first.
365  TiXmlElement *composite = field->FirstChildElement("C");
366 
367  while (composite)
368  {
369  nextCompositeNumber++;
370 
371  int indx;
372  int err = composite->QueryIntAttribute("ID", &indx);
373  ASSERTL0(err == TIXML_SUCCESS, "Unable to read attribute ID.");
374 // ASSERTL0(indx == nextCompositeNumber, "Composite IDs must begin with zero and be sequential.");
375 
376  TiXmlNode* compositeChild = composite->FirstChild();
377  // This is primarily to skip comments that may be present.
378  // Comments appear as nodes just like elements.
379  // We are specifically looking for text in the body
380  // of the definition.
381  while(compositeChild && compositeChild->Type() != TiXmlNode::TINYXML_TEXT)
382  {
383  compositeChild = compositeChild->NextSibling();
384  }
385 
386  ASSERTL0(compositeChild, "Unable to read composite definition body.");
387  std::string compositeStr = compositeChild->ToText()->ValueStr();
388 
389  /// Parse out the element components corresponding to type of element.
390 
391  std::istringstream compositeDataStrm(compositeStr.c_str());
392 
393  try
394  {
395  bool first = true;
396  std::string prevCompositeElementStr;
397 
398  while (!compositeDataStrm.fail())
399  {
400  std::string compositeElementStr;
401  compositeDataStrm >> compositeElementStr;
402 
403  if (!compositeDataStrm.fail())
404  {
405  if (first)
406  {
407  first = false;
408 
410  m_meshComposites[indx] = curVector;
411  }
412 
413  if (compositeElementStr.length() > 0)
414  {
415  ResolveGeomRef(prevCompositeElementStr, compositeElementStr, m_meshComposites[indx]);
416  }
417  prevCompositeElementStr = compositeElementStr;
418  }
419  }
420  }
421  catch(...)
422  {
424  (std::string("Unable to read COMPOSITE data for composite: ") + compositeStr).c_str());
425  }
426 
427  /// Keep looking
428  composite = composite->NextSiblingElement("C");
429  }
430  }
431 
432 
434  {
435  SegGeomSharedPtr returnval;
436  SegGeomMap::iterator x = m_segGeoms.find(eID);
437  ASSERTL0(x != m_segGeoms.end(), "Segment not found.");
438  return x->second;
439  };
440 
441 
442  // Take the string that is the composite reference and find the
443  // pointer to the Geometry object corresponding to it.
444 
445  // The only allowable combinations of previous and current items
446  // are V (0D); E (1D); and T and Q (2D). Only elements of the same
447  // dimension are allowed to be grouped.
448  void MeshGraph2D::ResolveGeomRef(const std::string &prevToken, const std::string &token,
449  Composite& composite)
450  {
451  try
452  {
453  std::istringstream tokenStream(token);
454  std::istringstream prevTokenStream(prevToken);
455 
456  char type;
457  char prevType;
458 
459  tokenStream >> type;
460 
461  std::string::size_type indxBeg = token.find_first_of('[') + 1;
462  std::string::size_type indxEnd = token.find_last_of(']') - 1;
463 
464  ASSERTL0(indxBeg <= indxEnd, (std::string("Error reading index definition:") + token).c_str());
465 
466  std::string indxStr = token.substr(indxBeg, indxEnd - indxBeg + 1);
467  std::vector<unsigned int> seqVector;
469 
470  bool err = ParseUtils::GenerateSeqVector(indxStr.c_str(), seqVector);
471 
472  ASSERTL0(err, (std::string("Error reading composite elements: ") + indxStr).c_str());
473 
474  prevTokenStream >> prevType;
475 
476  // All composites must be of the same dimension.
477  bool validSequence = (prevToken.empty() || // No previous, then current is just fine.
478  (type == 'V' && prevType == 'V') ||
479  (type == 'E' && prevType == 'E') ||
480  ((type == 'T' || type == 'Q') &&
481  (prevType == 'T' || prevType == 'Q')));
482 
483  ASSERTL0(validSequence, (std::string("Invalid combination of composite items: ")
484  + type + " and " + prevType + ".").c_str());
485 
486 
487  switch(type)
488  {
489  case 'E': // Edge
490  for (seqIter = seqVector.begin(); seqIter != seqVector.end(); ++seqIter)
491  {
492  if (m_segGeoms.find(*seqIter) == m_segGeoms.end())
493  {
494  char errStr[16] = "";
495  ::sprintf(errStr, "%d", *seqIter);
496  NEKERROR(ErrorUtil::ewarning, (std::string("Unknown edge index: ") + errStr).c_str());
497  }
498  else
499  {
500  composite->push_back(m_segGeoms[*seqIter]);
501  }
502  }
503  break;
504 
505  case 'T': // Triangle
506  {
507  for (seqIter = seqVector.begin(); seqIter != seqVector.end(); ++seqIter)
508  {
509  if (m_triGeoms.count(*seqIter) == 0 )
510  {
511  char errStr[16] = "";
512  ::sprintf(errStr, "%d", *seqIter);
513  NEKERROR(ErrorUtil::ewarning, (std::string("Unknown triangle index: ") + errStr).c_str());
514  }
515  else
516  {
517  if(CheckRange(*m_triGeoms[*seqIter]))
518  {
519  composite->push_back(m_triGeoms[*seqIter]);
520  }
521  }
522  }
523  }
524  break;
525 
526  case 'Q': // Quad
527  {
528  for (seqIter = seqVector.begin(); seqIter != seqVector.end(); ++seqIter)
529  {
530  if (m_quadGeoms.count(*seqIter) == 0)
531  {
532  char errStr[16] = "";
533  ::sprintf(errStr, "%d", *seqIter);
534  NEKERROR(ErrorUtil::ewarning, (std::string("Unknown quad index: ") + errStr +std::string(" in Composite section")).c_str());
535  }
536  else
537  {
538  if(CheckRange(*m_quadGeoms[*seqIter]))
539  {
540  composite->push_back(m_quadGeoms[*seqIter]);
541  }
542  }
543  }
544  }
545  break;
546 
547  case 'V': // Vertex
548  for (seqIter = seqVector.begin(); seqIter != seqVector.end(); ++seqIter)
549  {
550  if (*seqIter >= m_vertSet.size())
551  {
552  char errStr[16] = "";
553  ::sprintf(errStr, "%d", *seqIter);
554  NEKERROR(ErrorUtil::ewarning, (std::string("Unknown vertex index: ") + errStr).c_str());
555  }
556  else
557  {
558  composite->push_back(m_vertSet[*seqIter]);
559  }
560  }
561  break;
562 
563  default:
564  NEKERROR(ErrorUtil::efatal, (std::string("Unrecognized composite token: ") + token).c_str());
565  }
566  }
567  catch(...)
568  {
569  NEKERROR(ErrorUtil::efatal, (std::string("Problem processing composite token: ") + token).c_str());
570  }
571 
572  return;
573  }
574 
576  {
577  SegGeomSharedPtr Sedge;
578 
579  if(!(Sedge = boost::dynamic_pointer_cast<SegGeom>(edge)))
580  {
581  ASSERTL0(false,"Dynamics cast failed");
582 
583  }
584  return GetElementsFromEdge(Sedge);
585 
586  }
588  {
589  // Search tris and quads
590  // Need to iterate through vectors because there may be multiple
591  // occurrences.
592  ElementEdgeSharedPtr elementEdge;
593  //TriGeomVectorIter triIter;
594 
596 
597  CompositeMapIter compIter;
598  TriGeomSharedPtr triGeomShPtr;
599  QuadGeomSharedPtr quadGeomShPtr;
600 
601  GeometryVectorIter geomIter;
602 
603  for(int d = 0; d < m_domain.size(); ++d)
604  {
605  for (compIter = m_domain[d].begin(); compIter != m_domain[d].end(); ++compIter)
606  {
607  for (geomIter = (compIter->second)->begin(); geomIter != (compIter->second)->end(); ++geomIter)
608  {
609  triGeomShPtr = boost::dynamic_pointer_cast<TriGeom>(*geomIter);
610  quadGeomShPtr = boost::dynamic_pointer_cast<QuadGeom>(*geomIter);
611 
612  if (triGeomShPtr || quadGeomShPtr)
613  {
614  int edgeNum;
615  if (triGeomShPtr)
616  {
617  if ((edgeNum = triGeomShPtr->WhichEdge(edge)) > -1)
618  {
620  elementEdge->m_Element = triGeomShPtr;
621  elementEdge->m_EdgeIndx = edgeNum;
622  returnval->push_back(elementEdge);
623  }
624  }
625  else if (quadGeomShPtr)
626  {
627  if ((edgeNum = quadGeomShPtr->WhichEdge(edge)) > -1)
628  {
630  elementEdge->m_Element = quadGeomShPtr;
631  elementEdge->m_EdgeIndx = edgeNum;
632  returnval->push_back(elementEdge);
633  }
634  }
635  }
636  }
637  }
638  }
639 
640  return returnval;
641  }
642 
644  {
646  // Perhaps, a check should be done here to ensure that
647  // in case elements->size!=1, all elements to which
648  // the edge belongs have the same type and order of
649  // expansion such that no confusion can arise.
650  ExpansionShPtr expansion = GetExpansion((*elements)[0]->m_Element, variable);
651 
652  int edge_id = (*elements)[0]->m_EdgeIndx;
653 
654  if((*elements)[0]->m_Element->GetShapeType() == LibUtilities::eTriangle)
655  {
656  edge_id = (edge_id)? 1:0;
657  }
658  else
659  {
660  edge_id = edge_id%2;
661  }
662 
663  int nummodes = expansion->m_basisKeyVector[edge_id].GetNumModes();
664  int numpoints = expansion->m_basisKeyVector[edge_id].GetNumPoints();
665 
666  if((*elements)[0]->m_Element->GetShapeType() == LibUtilities::eTriangle)
667  {
668  // Use edge 0 to define basis of order relevant to edge
669  switch(expansion->m_basisKeyVector[edge_id].GetBasisType())
670  {
672  {
673  switch(expansion->m_basisKeyVector[edge_id].GetPointsType())
674  {
676  {
678  return LibUtilities::BasisKey(expansion->m_basisKeyVector[0].GetBasisType(),nummodes,pkey);
679  }
680  break;
681 
682  default:
683  ASSERTL0(false,"Unexpected points distribution");
684 
685  // It doesn't matter what we return
686  // here since the ASSERT will stop
687  // execution. Just return something
688  // to prevent warnings messages.
690  return LibUtilities::BasisKey(expansion->m_basisKeyVector[0].GetBasisType(),nummodes,pkey);
691  break;
692  }
693  }
694  break;
695  case LibUtilities::eOrtho_B: // Assume this is called from nodal triangular basis
696  {
697  switch(expansion->m_basisKeyVector[edge_id].GetPointsType())
698  {
700  {
703  }
704  break;
705 
706  default:
707  ASSERTL0(false,"Unexpected points distribution");
708 
709  // It doesn't matter what we return
710  // here since the ASSERT will stop
711  // execution. Just return something
712  // to prevent warnings messages.
714  return LibUtilities::BasisKey(expansion->m_basisKeyVector[0].GetBasisType(),nummodes,pkey);
715  break;
716  }
717  }
718  break;
720  {
721  switch(expansion->m_basisKeyVector[edge_id].GetPointsType())
722  {
724  {
726  return LibUtilities::BasisKey(expansion->m_basisKeyVector[0].GetBasisType(),nummodes,pkey);
727  }
728  break;
730  {
732  return LibUtilities::BasisKey(expansion->m_basisKeyVector[0].GetBasisType(),nummodes,pkey);
733  }
734  break;
735 
736  default:
737  ASSERTL0(false,"Unexpected points distribution");
738 
739  // It doesn't matter what we return
740  // here since the ASSERT will stop
741  // execution. Just return something
742  // to prevent warnings messages.
744  return LibUtilities::BasisKey(expansion->m_basisKeyVector[0].GetBasisType(),nummodes,pkey);
745  break;
746  }
747  }
748  break;
750  {
751  switch(expansion->m_basisKeyVector[edge_id].GetPointsType())
752  {
754  {
756  return LibUtilities::BasisKey(expansion->m_basisKeyVector[0].GetBasisType(),nummodes,pkey);
757  }
758  break;
759  default:
760  ASSERTL0(false,"Unexpected points distribution");
761  // It doesn't matter what we return here
762  // since the ASSERT will stop execution.
763  // Just return something to prevent
764  // warnings messages.
766  return LibUtilities::BasisKey(expansion->m_basisKeyVector[0].GetBasisType(),nummodes,pkey);
767  break;
768  }
769  }
770  break;
771 
772  default:
773  ASSERTL0(false,"Unexpected basis distribution");
774  // It doesn't matter what we return here since the
775  // ASSERT will stop execution. Just return
776  // something to prevent warnings messages.
778  return LibUtilities::BasisKey(expansion->m_basisKeyVector[0].GetBasisType(),nummodes,pkey);
779  }
780  }
781  else
782  {
783  // Quadrilateral
784  const LibUtilities::PointsKey pkey(numpoints,expansion->m_basisKeyVector[edge_id].GetPointsType());
785  return LibUtilities::BasisKey(expansion->m_basisKeyVector[edge_id].GetBasisType(),nummodes,pkey);
786  }
787 
788  ASSERTL0(false, "Unable to determine edge points type.");
790  }
791  }; //end of namespace
792 }; //end of namespace
793 
794 //
795 // $Log: MeshGraph2D.cpp,v $
796 // Revision 1.41 2009/12/17 02:08:04 bnelson
797 // Fixed visual studio compiler warning.
798 //
799 // Revision 1.40 2009/10/22 17:41:47 cbiotto
800 // Update for variable order expansion
801 //
802 // Revision 1.39 2009/07/02 13:32:24 sehunchun
803 // *** empty log message ***
804 //
805 // Revision 1.38 2009/07/02 11:39:44 sehunchun
806 // Modification for 2D gemoetry embedded in 3D
807 //
808 // Revision 1.37 2009/05/01 13:23:21 pvos
809 // Fixed various bugs
810 //
811 // Revision 1.36 2009/04/20 16:13:23 sherwin
812 // Modified Import and Write functions and redefined how Expansion is used
813 //
814 // Revision 1.35 2009/01/21 16:59:03 pvos
815 // Added additional geometric factors to improve efficiency
816 //
817 // Revision 1.34 2009/01/12 10:26:59 pvos
818 // Added input tags for nodal expansions
819 //
820 // Revision 1.33 2008/09/12 11:26:19 pvos
821 // Updates for mappings in 3D
822 //
823 // Revision 1.32 2008/09/09 14:21:44 sherwin
824 // Updates for first working version of curved edges
825 //
826 // Revision 1.31 2008/08/14 22:11:03 sherwin
827 // Mods for HDG update
828 //
829 // Revision 1.30 2008/07/29 22:23:36 sherwin
830 // various mods for DG advection solver in Multiregions. Added virtual calls to Geometry, Geometry1D, 2D and 3D
831 //
832 // Revision 1.29 2008/05/29 19:07:39 delisi
833 // Removed the Write(..) methods, so it is only in the base MeshGraph class. Also, added a line to set the global ID of the geometry object for every element read in.
834 //
835 // Revision 1.28 2008/05/28 21:42:18 jfrazier
836 // Minor comment spelling change.
837 //
838 // Revision 1.27 2008/03/18 14:14:49 pvos
839 // Update for nodal triangular helmholtz solver
840 //
841 // Revision 1.26 2008/03/11 17:02:24 jfrazier
842 // Modified the GetElementsFromEdge to use the domain.
843 //
844 // Revision 1.25 2008/01/21 19:58:14 sherwin
845 // Updated so that QuadGeom and TriGeom have SegGeoms instead of EdgeComponents
846 //
847 // Revision 1.24 2007/12/17 20:27:24 sherwin
848 // Added normals to GeomFactors
849 //
850 // Revision 1.23 2007/12/14 17:39:18 jfrazier
851 // Fixed composite reader to handle ranges and comma separated lists.
852 //
853 // Revision 1.22 2007/12/11 21:51:53 jfrazier
854 // Updated 2d components so elements could be retrieved from edges.
855 //
856 // Revision 1.21 2007/12/04 03:02:26 jfrazier
857 // Changed to stringstream.
858 //
859 // Revision 1.20 2007/09/20 22:25:06 jfrazier
860 // Added expansion information to meshgraph class.
861 //
862 // Revision 1.19 2007/07/26 01:38:33 jfrazier
863 // Cleanup of some attribute reading code.
864 //
865 // Revision 1.18 2007/07/24 16:52:09 jfrazier
866 // Added domain code.
867 //
868 // Revision 1.17 2007/07/05 04:21:10 jfrazier
869 // Changed id format and propagated from 1d to 2d.
870 //
871 // Revision 1.16 2007/06/10 02:27:11 jfrazier
872 // Another checkin with an incremental completion of the boundary conditions reader
873 //
874 // Revision 1.15 2007/06/07 23:55:24 jfrazier
875 // Intermediate revisions to add parsing for boundary conditions file.
876 //
877 // Revision 1.14 2007/05/28 21:48:42 sherwin
878 // Update for 2D functionality
879 //
880 // Revision 1.13 2006/10/17 22:26:01 jfrazier
881 // Added capability to specify ranges in composite definition.
882 //
883 // Revision 1.12 2006/10/17 18:42:54 jfrazier
884 // Removed "NUMBER" attribute in items.
885 //
886 // Revision 1.11 2006/09/26 23:41:53 jfrazier
887 // Updated to account for highest level NEKTAR tag and changed the geometry tag to GEOMETRY.
888 //
889 // Revision 1.10 2006/08/24 18:50:00 jfrazier
890 // Completed error checking on permissable composite item combinations.
891 //
892 // Revision 1.9 2006/08/18 19:37:17 jfrazier
893 // *** empty log message ***
894 //
895 // Revision 1.8 2006/08/17 22:55:00 jfrazier
896 // Continued adding code to process composites in the mesh2d.
897 //
898 // Revision 1.7 2006/08/16 23:34:42 jfrazier
899 // *** empty log message ***
900 //
901 // Revision 1.6 2006/06/02 18:48:40 sherwin
902 // Modifications to make ProjectLoc2D run bit there are bus errors for order > 3
903 //
904 // Revision 1.5 2006/06/01 14:15:30 sherwin
905 // Added typdef of boost wrappers and made GeoFac a boost shared pointer.
906 //
907 // Revision 1.4 2006/05/30 14:00:04 sherwin
908 // Updates to make MultiRegions and its Demos work
909 //
910 // Revision 1.3 2006/05/23 19:56:33 jfrazier
911 // These build and run, but the expansion pieces are commented out
912 // because they would not run.
913 //
914 // Revision 1.2 2006/05/09 13:37:01 jfrazier
915 // Removed duplicate definition of shared vertex pointer.
916 //
917 // Revision 1.1 2006/05/04 18:59:01 kirby
918 // *** empty log message ***
919 //
920 // Revision 1.11 2006/04/11 23:18:11 jfrazier
921 // Completed MeshGraph2D for tri's and quads. Not thoroughly tested.
922 //
923 // Revision 1.10 2006/04/09 02:08:35 jfrazier
924 // Added precompiled header.
925 //
926 // Revision 1.9 2006/04/04 23:12:37 jfrazier
927 // More updates to readers. Still work to do on MeshGraph2D to store tris and quads.
928 //
929 // Revision 1.8 2006/03/25 00:58:29 jfrazier
930 // Many changes dealing with fundamental structure and reading/writing.
931 //
932 // Revision 1.7 2006/03/12 14:20:43 sherwin
933 //
934 // First compiling version of SpatialDomains and associated modifications
935 //
936 // Revision 1.6 2006/03/12 07:42:03 sherwin
937 //
938 // Updated member names and StdRegions call. Still has not been compiled
939 //
940 // Revision 1.5 2006/02/26 21:19:43 bnelson
941 // Fixed a variety of compiler errors caused by updates to the coding standard.
942 //
943 // Revision 1.4 2006/02/19 01:37:34 jfrazier
944 // Initial attempt at bringing into conformance with the coding standard. Still more work to be done. Has not been compiled.
945 //
946 //
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:135
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:132
void ReadEdges(TiXmlDocument &doc)
Definition: MeshGraph2D.cpp:97
void ReadComposites(TiXmlDocument &doc)
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
void ReadCurves(TiXmlDocument &doc)
Definition: MeshGraph.cpp:1103
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
boost::shared_ptr< QuadGeom > QuadGeomSharedPtr
Definition: QuadGeom.h:54
void ReadElements(TiXmlDocument &doc)
LibUtilities::BasisKey GetEdgeBasisKey(SegGeomSharedPtr edge, const std::string variable="DefaultVar")
Return the BasisKey corresponding to an edge of an element If the expansion is a triangle the Modifie...
Principle Modified Functions .
Definition: BasisType.h:49
PointGeomSharedPtr GetVertex(int id)
Definition: MeshGraph.h:567
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:50
static bool GenerateSeqVector(const char *const str, std::vector< unsigned int > &vec)
Definition: ParseUtils.hpp:78
virtual void ReadGeometry(const std::string &infilename)
Read will read the meshgraph vertices given a filename.
Definition: MeshGraph.cpp:230
static StdRegions::Orientation GetEdgeOrientation(const SegGeom &edge1, const SegGeom &edge2)
Get the orientation of edge1.
Definition: SegGeom.cpp:293
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:57
std::vector< GeometrySharedPtr >::iterator GeometryVectorIter
Definition: Geometry.h:58
Principle Orthogonal Functions .
Definition: BasisType.h:47
boost::shared_ptr< ElementEdge > ElementEdgeSharedPtr
Definition: MeshGraph.h:129
void ReadGeometry(const std::string &infilename)
Read will read the meshgraph vertices given a filename.
Definition: MeshGraph2D.cpp:62
boost::shared_ptr< SegGeom > SegGeomSharedPtr
Definition: Geometry2D.h:60
Principle Modified Functions .
Definition: BasisType.h:50
boost::shared_ptr< DomainRange > DomainRangeShPtr
Definition: MeshGraph.h:154
std::map< int, Composite >::iterator CompositeMapIter
Definition: MeshGraph.h:113
Defines a specification for a set of points.
Definition: Points.h:58
boost::shared_ptr< GeometryVector > Composite
Definition: MeshGraph.h:111
ElementEdgeVectorSharedPtr GetElementsFromEdge(SegGeomSharedPtr edge)
Return the elements (shared ptrs) that have this edge.
static const int kNedges
Get the orientation of face1.
Definition: TriGeom.h:98
boost::shared_ptr< ElementEdgeVector > ElementEdgeVectorSharedPtr
Definition: MeshGraph.h:131
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
std::vector< CompositeMap > m_domain
Definition: MeshGraph.h:422
boost::shared_ptr< Expansion > ExpansionShPtr
Definition: MeshGraph.h:170
Base class for a spectral/hp element mesh.
Definition: MeshGraph.h:183
void ReadDomain(TiXmlDocument &doc)
Definition: MeshGraph.cpp:1009
bool CheckRange(Geometry2D &geom)
Check if goemetry is in range definition if activated.
Definition: MeshGraph.cpp:1813
boost::shared_ptr< Geometry1D > Geometry1DSharedPtr
Definition: Geometry1D.h:48
ExpansionShPtr GetExpansion(GeometrySharedPtr geom, const std::string variable="DefaultVar")
Definition: MeshGraph.cpp:2136
boost::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:58
static const BasisKey NullBasisKey(eNoBasisType, 0, NullPointsKey)
Defines a null basis with no type or points.
void ReadExpansions(const std::string &infilename)
Read the expansions given the XML file path.
Definition: MeshGraph.cpp:525
Lagrange for SEM basis .
Definition: BasisType.h:53
void ResolveGeomRef(const std::string &prevToken, const std::string &token, Composite &composite)
Describes the specification for a Basis.
Definition: Basis.h:50
boost::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:60
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50
SegGeomSharedPtr GetSegGeom(int eID)