Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | List of all members
Nektar::SpatialDomains::MeshGraph3D Class Reference

#include <MeshGraph3D.h>

Inheritance diagram for Nektar::SpatialDomains::MeshGraph3D:
Inheritance graph
[legend]
Collaboration diagram for Nektar::SpatialDomains::MeshGraph3D:
Collaboration graph
[legend]

Public Member Functions

 MeshGraph3D ()
 MeshGraph3D (const LibUtilities::SessionReaderSharedPtr &pSession, const DomainRangeShPtr &rng=NullDomainRangeShPtr)
virtual ~MeshGraph3D ()
void ReadGeometry (const std::string &infilename)
 Read will read the meshgraph vertices given a filename.
void ReadGeometry (TiXmlDocument &doc)
 Read will read the meshgraph vertices given a TiXmlDocument.
SegGeomSharedPtr GetSegGeom (int eID)
Geometry2DSharedPtr GetGeometry2D (int gID)
int GetCoordim (void)
const TriGeomMapGetTrigeoms (void) const
const QuadGeomMapGetQuadgeoms (void) const
void GenXGeoFac ()
int GetNseggeoms () const
int GetVidFromElmt (LibUtilities::ShapeType shape, const int vert, const int elmt) const
int GetEidFromElmt (LibUtilities::ShapeType shape, const int edge, const int elmt) const
StdRegions::Orientation GetEorientFromElmt (LibUtilities::ShapeType shape, const int edge, const int elmt) const
StdRegions::Orientation GetCartesianEorientFromElmt (LibUtilities::ShapeType shape, const int edge, const int elmt) const
int GetNumComposites (void)
int GetNumCompositeItems (int whichComposite)
ElementFaceVectorSharedPtr GetElementsFromFace (Geometry2DSharedPtr face)
 Return the elements (shared ptrs) that have this face.
LibUtilities::BasisKey GetFaceBasisKey (Geometry2DSharedPtr face, const int flag, const std::string variable="DefaultVar")
 Return the BasisKey corresponding to a face of an element.
- Public Member Functions inherited from Nektar::SpatialDomains::MeshGraph
 MeshGraph ()
 MeshGraph (unsigned int meshDimension, unsigned int spaceDimension)
 MeshGraph (const LibUtilities::SessionReaderSharedPtr &pSession, const DomainRangeShPtr &rng=NullDomainRangeShPtr)
virtual ~MeshGraph ()
void ReadGeometryInfo (const std::string &infilename)
 Read geometric information from a file.
void ReadGeometryInfo (TiXmlDocument &doc)
 Read geometric information from an XML document.
void ReadExpansions (const std::string &infilename)
 Read the expansions given the XML file path.
void ReadExpansions (TiXmlDocument &doc)
 Read the expansions given the XML document reference.
void ReadDomain (TiXmlDocument &doc)
void ReadCurves (TiXmlDocument &doc)
void ReadCurves (std::string &infilename)
int GetMeshDimension () const
 Dimension of the mesh (can be a 1D curve in 3D space).
int GetSpaceDimension () const
 Dimension of the space (can be a 1D curve in 3D space).
void SetDomainRange (NekDouble xmin, NekDouble xmax, NekDouble ymin=NekConstants::kNekUnsetDouble, NekDouble ymax=NekConstants::kNekUnsetDouble, NekDouble zmin=NekConstants::kNekUnsetDouble, NekDouble zmax=NekConstants::kNekUnsetDouble)
bool CheckRange (Geometry2D &geom)
 Check if goemetry is in range definition if activated.
bool CheckRange (Geometry3D &geom)
 Check if goemetry is in range definition if activated.
Composite GetComposite (int whichComposite) const
GeometrySharedPtr GetCompositeItem (int whichComposite, int whichItem)
void GetCompositeList (const std::string &compositeStr, CompositeMap &compositeVector) const
const CompositeMapGetComposites () const
const std::vector< CompositeMap > & GetDomain (void) const
const CompositeMapGetDomain (int domain) const
const ExpansionMapGetExpansions ()
const ExpansionMapGetExpansions (const std::string variable)
ExpansionShPtr GetExpansion (GeometrySharedPtr geom, const std::string variable="DefaultVar")
void SetExpansions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
 Sets expansions given field definitions.
void SetExpansions (std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef, std::vector< std::vector< LibUtilities::PointsType > > &pointstype)
 Sets expansions given field definition, quadrature points.
void SetExpansionsToEvenlySpacedPoints (int npoints=0)
 Sets expansions to have equispaced points.
void SetExpansions (const std::string variable, ExpansionMapShPtr &exp)
 This function sets the expansion #exp in map with entry #variable.
void SetBasisKey (LibUtilities::ShapeType shape, LibUtilities::BasisKeyVector &keys, std::string var="DefaultVar")
 Sets the basis key for all expansions of the given shape.
bool SameExpansions (const std::string var1, const std::string var2)
bool CheckForGeomInfo (std::string parameter)
const std::string GetGeomInfo (std::string parameter)
LibUtilities::BasisKeyVector DefineBasisKeyFromExpansionTypeHomo (GeometrySharedPtr in, ExpansionType type_x, ExpansionType type_y, ExpansionType type_z, const int nummodes_x, const int nummodes_y, const int nummodes_z)
int GetNvertices () const
PointGeomSharedPtr GetVertex (int id)
PointGeomSharedPtr AddVertex (NekDouble x, NekDouble y, NekDouble z)
 Adds a vertex to the with the next available ID.
SegGeomSharedPtr AddEdge (PointGeomSharedPtr v0, PointGeomSharedPtr v1, CurveSharedPtr curveDefinition=CurveSharedPtr())
 Adds an edge between two points. If curveDefinition is null, then the edge is straight, otherwise it is curved according to the curveDefinition.
SegGeomSharedPtr GetEdge (unsigned int id)
TriGeomSharedPtr AddTriangle (SegGeomSharedPtr edges[], StdRegions::Orientation orient[])
QuadGeomSharedPtr AddQuadrilateral (SegGeomSharedPtr edges[], StdRegions::Orientation orient[])
TetGeomSharedPtr AddTetrahedron (TriGeomSharedPtr tfaces[TetGeom::kNtfaces])
PyrGeomSharedPtr AddPyramid (TriGeomSharedPtr tfaces[PyrGeom::kNtfaces], QuadGeomSharedPtr qfaces[PyrGeom::kNqfaces])
PrismGeomSharedPtr AddPrism (TriGeomSharedPtr tfaces[PrismGeom::kNtfaces], QuadGeomSharedPtr qfaces[PrismGeom::kNqfaces])
HexGeomSharedPtr AddHexahedron (QuadGeomSharedPtr qfaces[HexGeom::kNqfaces])
const CurveVectorGetCurvedEdges () const
const CurveVectorGetCurvedFaces () const
const SegGeomMapGetAllSegGeoms () const
const TriGeomMapGetAllTriGeoms () const
const QuadGeomMapGetAllQuadGeoms () const
const TetGeomMapGetAllTetGeoms () const
const PyrGeomMapGetAllPyrGeoms () const
const PrismGeomMapGetAllPrismGeoms () const
const HexGeomMapGetAllHexGeoms () const
template<typename ElementType >
const std::map< int,
boost::shared_ptr< ElementType > > & 
GetAllElementsOfType () const
 Convenience method for ElVis.

Protected Member Functions

void ReadEdges (TiXmlDocument &doc)
void ReadFaces (TiXmlDocument &doc)
void ReadElements (TiXmlDocument &doc)
void ReadComposites (TiXmlDocument &doc)
void ResolveGeomRef (const std::string &prevToken, const std::string &token, Composite &composite)
- Protected Member Functions inherited from Nektar::SpatialDomains::MeshGraph
ExpansionMapShPtr SetUpExpansionMap (void)

Private Member Functions

void PopulateFaceToElMap (Geometry3DSharedPtr element, int kNfaces)
 Given a 3D geometry object #element, populate the face to element map m_faceToElMap which maps faces to their corresponding element(s).

Private Attributes

boost::unordered_map< int,
ElementFaceVectorSharedPtr
m_faceToElMap

Additional Inherited Members

- Static Public Member Functions inherited from Nektar::SpatialDomains::MeshGraph
static boost::shared_ptr
< MeshGraph
Read (const LibUtilities::SessionReaderSharedPtr &pSession, DomainRangeShPtr &rng=NullDomainRangeShPtr)
static boost::shared_ptr
< MeshGraph
Read (const std::string &infilename, bool pReadExpansions=true)
static LibUtilities::BasisKeyVector DefineBasisKeyFromExpansionType (GeometrySharedPtr in, ExpansionType type, const int order)
- Protected Attributes inherited from Nektar::SpatialDomains::MeshGraph
LibUtilities::SessionReaderSharedPtr m_session
PointGeomMap m_vertSet
InterfaceCompList m_iComps
CurveVector m_curvedEdges
CurveVector m_curvedFaces
SegGeomMap m_segGeoms
TriGeomMap m_triGeoms
QuadGeomMap m_quadGeoms
TetGeomMap m_tetGeoms
PyrGeomMap m_pyrGeoms
PrismGeomMap m_prismGeoms
HexGeomMap m_hexGeoms
int m_meshDimension
int m_spaceDimension
int m_partition
bool m_meshPartitioned
CompositeMap m_meshComposites
std::vector< CompositeMapm_domain
DomainRangeShPtr m_domainRange
ExpansionMapShPtrMap m_expansionMapShPtrMap
GeomInfoMap m_geomInfo

Detailed Description

Definition at line 55 of file MeshGraph3D.h.

Constructor & Destructor Documentation

Nektar::SpatialDomains::MeshGraph3D::MeshGraph3D ( )

Definition at line 46 of file MeshGraph3D.cpp.

: MeshGraph(3,3)
{
}
Nektar::SpatialDomains::MeshGraph3D::MeshGraph3D ( const LibUtilities::SessionReaderSharedPtr pSession,
const DomainRangeShPtr rng = NullDomainRangeShPtr 
)

Definition at line 50 of file MeshGraph3D.cpp.

References Nektar::SpatialDomains::MeshGraph::ReadExpansions(), and ReadGeometry().

: MeshGraph(pSession,rng)
{
ReadGeometry(pSession->GetDocument());
ReadExpansions(pSession->GetDocument());
}
Nektar::SpatialDomains::MeshGraph3D::~MeshGraph3D ( )
virtual

Definition at line 58 of file MeshGraph3D.cpp.

{
}

Member Function Documentation

void Nektar::SpatialDomains::MeshGraph3D::GenXGeoFac ( )
StdRegions::Orientation Nektar::SpatialDomains::MeshGraph3D::GetCartesianEorientFromElmt ( LibUtilities::ShapeType  shape,
const int  edge,
const int  elmt 
) const
inline

Definition at line 147 of file MeshGraph3D.h.

References ASSERTL2, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eForwards, Nektar::LibUtilities::eTriangle, Nektar::SpatialDomains::MeshGraph::m_quadGeoms, and Nektar::SpatialDomains::MeshGraph::m_triGeoms.

{
{
ASSERTL2(m_triGeoms.find(elmt) != m_triGeoms.end(),
"eid is out of range");
returnval = m_triGeoms.find(elmt)->second->GetEorient(edge);
}
else
{
ASSERTL2(m_quadGeoms.find(elmt) != m_quadGeoms.end(),
"eid is out of range");
returnval = m_quadGeoms.find(elmt)->second->GetEorient(edge);
}
// swap orientation if on edge 2 & 3 (if quad)
if(edge >= 2)
{
if(returnval == StdRegions::eForwards)
{
}
else
{
returnval = StdRegions::eForwards;
}
}
return returnval;
}
int Nektar::SpatialDomains::MeshGraph3D::GetCoordim ( void  )
inline

Definition at line 69 of file MeshGraph3D.h.

References Nektar::SpatialDomains::MeshGraph::GetSpaceDimension().

{
}
int Nektar::SpatialDomains::MeshGraph3D::GetEidFromElmt ( LibUtilities::ShapeType  shape,
const int  edge,
const int  elmt 
) const
inline

Definition at line 109 of file MeshGraph3D.h.

References ASSERTL2, Nektar::LibUtilities::eTriangle, Nektar::SpatialDomains::MeshGraph::m_quadGeoms, and Nektar::SpatialDomains::MeshGraph::m_triGeoms.

{
{
ASSERTL2(m_triGeoms.find(elmt) != m_triGeoms.end(),
"eid is out of range");
return m_triGeoms.find(elmt)->second->GetEid(edge);
}
else
{
ASSERTL2(m_quadGeoms.find(elmt) != m_quadGeoms.end(),
"eid is out of range");
return m_quadGeoms.find(elmt)->second->GetEid(edge);
}
}
ElementFaceVectorSharedPtr Nektar::SpatialDomains::MeshGraph3D::GetElementsFromFace ( Geometry2DSharedPtr  face)

Return the elements (shared ptrs) that have this face.

Definition at line 986 of file MeshGraph3D.cpp.

References ASSERTL0, Nektar::iterator, and m_faceToElMap.

Referenced by GetFaceBasisKey().

{
m_faceToElMap.find(face->GetGlobalID());
ASSERTL0(it != m_faceToElMap.end(), "Unable to find corresponding face!");
return it->second;
}
StdRegions::Orientation Nektar::SpatialDomains::MeshGraph3D::GetEorientFromElmt ( LibUtilities::ShapeType  shape,
const int  edge,
const int  elmt 
) const
inline

Definition at line 128 of file MeshGraph3D.h.

References ASSERTL2, Nektar::LibUtilities::eTriangle, Nektar::SpatialDomains::MeshGraph::m_quadGeoms, and Nektar::SpatialDomains::MeshGraph::m_triGeoms.

{
{
ASSERTL2(m_triGeoms.find(elmt) != m_triGeoms.end(),
"eid is out of range");
return m_triGeoms.find(elmt)->second->GetEorient(edge);
}
else
{
ASSERTL2(m_quadGeoms.find(elmt) != m_quadGeoms.end(),
"eid is out of range");
return m_quadGeoms.find(elmt)->second->GetEorient(edge);
}
}
LibUtilities::BasisKey Nektar::SpatialDomains::MeshGraph3D::GetFaceBasisKey ( Geometry2DSharedPtr  face,
const int  facedir,
const std::string  variable = "DefaultVar" 
)

Return the BasisKey corresponding to a face of an element.

Retrieve the basis key for a given face direction.

Definition at line 999 of file MeshGraph3D.cpp.

References ASSERTL0, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eGLL_Lagrange, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eModified_C, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::SpatialDomains::Geometry3D::GetDir(), GetElementsFromFace(), Nektar::SpatialDomains::MeshGraph::GetExpansion(), and Nektar::LibUtilities::NullBasisKey().

{
// Retrieve the list of elements and the associated face index
// to which the face geometry belongs.
ASSERTL0(elements->size() > 0, "No elements for the given face."
" Check all elements belong to the domain composite.");
// Perhaps, a check should be done here to ensure that in case
// elements->size!=1, all elements to which the edge belongs have
// the same type and order of expansion such that no confusion can
// arise.
// Get the Expansion structure detailing the basis keys used for
// this element.
ExpansionShPtr expansion = GetExpansion((*elements)[0]->m_Element,
variable);
// Retrieve the geometry object of the element as a Geometry3D.
boost::dynamic_pointer_cast<SpatialDomains::Geometry3D>(
expansion->m_geomShPtr);
// Use the geometry of the element to calculate the coordinate
// direction of the element which corresponds to the requested
// coordinate direction of the given face.
int dir = geom3d->GetDir((*elements)[0]->m_FaceIndx, facedir);
// Obtain the number of modes for the element basis key in this
// direction.
int nummodes = (int) expansion->m_basisKeyVector[dir].GetNumModes();
int numpoints = (int) expansion->m_basisKeyVector[dir].GetNumPoints();
switch(expansion->m_basisKeyVector[dir].GetBasisType())
{
{
switch (facedir)
{
case 0:
{
const LibUtilities::PointsKey pkey(nummodes+1,LibUtilities::eGaussLobattoLegendre);
return LibUtilities::BasisKey(LibUtilities::eModified_A,nummodes,pkey);
}
break;
case 1:
{
const LibUtilities::PointsKey pkey(nummodes+1,LibUtilities::eGaussLobattoLegendre);
if (face->GetNumVerts() == 3)
{
// Triangle
return LibUtilities::BasisKey(LibUtilities::eModified_B,nummodes,pkey);
}
else {
// Quadrilateral
return LibUtilities::BasisKey(LibUtilities::eModified_A,nummodes,pkey);
}
}
break;
default:
ASSERTL0(false,"invalid value to flag");
break;
}
}
break;
{
TriGeomSharedPtr triangle = boost::dynamic_pointer_cast<TriGeom>(face);
QuadGeomSharedPtr quadrilateral = boost::dynamic_pointer_cast<QuadGeom>(face);
if(quadrilateral)
{
const LibUtilities::PointsKey pkey(numpoints,LibUtilities::eGaussLobattoLegendre);
return LibUtilities::BasisKey(LibUtilities::eGLL_Lagrange,nummodes,pkey);
}
else if(triangle)
{
switch (facedir)
{
case 0:
{
const LibUtilities::PointsKey pkey(nummodes+1,LibUtilities::eGaussLobattoLegendre);
return LibUtilities::BasisKey(LibUtilities::eOrtho_A,nummodes,pkey);
}
break;
case 1:
{
const LibUtilities::PointsKey pkey(nummodes,LibUtilities::eGaussRadauMAlpha1Beta0);
return LibUtilities::BasisKey(LibUtilities::eOrtho_B,nummodes,pkey);
}
break;
default:
ASSERTL0(false,"invalid value to flag");
break;
}
}
else
{
ASSERTL0(false,"dynamic cast to a proper Geometry2D failed");
}
}
break;
{
switch (facedir)
{
case 0:
{
const LibUtilities::PointsKey pkey(nummodes+1,LibUtilities::eGaussLobattoLegendre);
return LibUtilities::BasisKey(LibUtilities::eOrtho_A,nummodes,pkey);
}
break;
case 1:
{
const LibUtilities::PointsKey pkey(nummodes,LibUtilities::eGaussRadauMAlpha1Beta0);
return LibUtilities::BasisKey(LibUtilities::eOrtho_B,nummodes,pkey);
}
break;
default:
ASSERTL0(false,"invalid value to flag");
break;
}
}
break;
// case eGLL_Lagrange_SEM:
// {
// const LibUtilities::PointsKey pkey(nummodes,LibUtilities::eGaussLobattoLegendre);
// return LibUtilities::BasisKey(LibUtilities::eGLL_Lagrange,nummodes,pkey);
// }
// break;
default:
ASSERTL0(false,"expansion type unknown");
break;
}
return LibUtilities::NullBasisKey; // Keep things happy by returning a value.
}
Geometry2DSharedPtr Nektar::SpatialDomains::MeshGraph3D::GetGeometry2D ( int  gID)

Definition at line 733 of file MeshGraph3D.cpp.

References Nektar::SpatialDomains::MeshGraph::m_quadGeoms, and Nektar::SpatialDomains::MeshGraph::m_triGeoms.

Referenced by ReadElements(), and ResolveGeomRef().

{
it1 = m_triGeoms.find(gID);
if (it1 != m_triGeoms.end())
return it1->second;
it2 = m_quadGeoms.find(gID);
if (it2 != m_quadGeoms.end())
return it2->second;
};
int Nektar::SpatialDomains::MeshGraph3D::GetNseggeoms ( ) const
inline

Definition at line 85 of file MeshGraph3D.h.

References Nektar::SpatialDomains::MeshGraph::m_segGeoms.

{
return int(m_segGeoms.size());
}
int Nektar::SpatialDomains::MeshGraph3D::GetNumCompositeItems ( int  whichComposite)
inline

Definition at line 186 of file MeshGraph3D.h.

References ErrorUtil::efatal, Nektar::SpatialDomains::MeshGraph::m_meshComposites, and NEKERROR.

{
int returnval = -1;
try
{
returnval = int(m_meshComposites[whichComposite]->size());
}
catch(...)
{
std::ostringstream errStream;
errStream << "Unable to access composite item [" << whichComposite << "].";
NEKERROR(ErrorUtil::efatal, errStream.str());
}
return returnval;
}
int Nektar::SpatialDomains::MeshGraph3D::GetNumComposites ( void  )
inline

Definition at line 181 of file MeshGraph3D.h.

References Nektar::SpatialDomains::MeshGraph::m_meshComposites.

{
return int(m_meshComposites.size());
}
const QuadGeomMap& Nektar::SpatialDomains::MeshGraph3D::GetQuadgeoms ( void  ) const
inline

Definition at line 78 of file MeshGraph3D.h.

References Nektar::SpatialDomains::MeshGraph::m_quadGeoms.

{
return m_quadGeoms;
}
SegGeomSharedPtr Nektar::SpatialDomains::MeshGraph3D::GetSegGeom ( int  eID)

Definition at line 724 of file MeshGraph3D.cpp.

References ASSERTL0, Nektar::iterator, and Nektar::SpatialDomains::MeshGraph::m_segGeoms.

Referenced by ReadFaces().

{
SegGeomSharedPtr returnval;
ASSERTL0(x != m_segGeoms.end(), "Segment "
+ boost::lexical_cast<string>(eID) + " not found.");
return x->second;
};
const TriGeomMap& Nektar::SpatialDomains::MeshGraph3D::GetTrigeoms ( void  ) const
inline

Definition at line 73 of file MeshGraph3D.h.

References Nektar::SpatialDomains::MeshGraph::m_triGeoms.

{
return m_triGeoms;
}
int Nektar::SpatialDomains::MeshGraph3D::GetVidFromElmt ( LibUtilities::ShapeType  shape,
const int  vert,
const int  elmt 
) const
inline

Definition at line 90 of file MeshGraph3D.h.

References ASSERTL2, Nektar::LibUtilities::eTriangle, Nektar::SpatialDomains::MeshGraph::m_quadGeoms, and Nektar::SpatialDomains::MeshGraph::m_triGeoms.

{
{
ASSERTL2(m_triGeoms.find(elmt) != m_triGeoms.end(),
"eid is out of range");
return m_triGeoms.find(elmt)->second->GetVid(vert);
}
else
{
ASSERTL2(m_quadGeoms.find(elmt) != m_quadGeoms.end(),
"eid is out of range");
return m_quadGeoms.find(elmt)->second->GetVid(vert);
}
}
void Nektar::SpatialDomains::MeshGraph3D::PopulateFaceToElMap ( Geometry3DSharedPtr  element,
int  kNfaces 
)
private

Given a 3D geometry object #element, populate the face to element map m_faceToElMap which maps faces to their corresponding element(s).

Parameters
elementElement to process.
kNfacesNumber of faces of #element. Should be removed and put into Geometry3D as a virtual member function.

Definition at line 1152 of file MeshGraph3D.cpp.

References Nektar::iterator, and m_faceToElMap.

Referenced by ReadElements().

{
// Set up face -> element map
for (int i = 0; i < kNfaces; ++i)
{
int faceId = element->GetFace(i)->GetGlobalID();
ElementFaceSharedPtr elementFace =
MemoryManager<ElementFace>::AllocateSharedPtr();
elementFace->m_Element = element;
elementFace->m_FaceIndx = i;
// Search map to see if face already exists.
m_faceToElMap.find(faceId);
if (it == m_faceToElMap.end())
{
MemoryManager<ElementFaceVector>::AllocateSharedPtr();
tmp->push_back(elementFace);
m_faceToElMap[faceId] = tmp;
}
else
{
ElementFaceVectorSharedPtr tmp = it->second;
tmp->push_back(elementFace);
}
}
}
void Nektar::SpatialDomains::MeshGraph3D::ReadComposites ( TiXmlDocument &  doc)
protected

We know we have it since we made it this far.

Look for elements in ELEMENT block.

All elements are of the form: "<C ID = "N"> ... </C>".

Read the ID field first.

Parse out the element components corresponding to type of element.

Keep looking

Definition at line 636 of file MeshGraph3D.cpp.

References ASSERTL0, ErrorUtil::efatal, Nektar::SpatialDomains::MeshGraph::m_meshComposites, NEKERROR, and ResolveGeomRef().

Referenced by ReadGeometry().

{
TiXmlHandle docHandle(&doc);
/// We know we have it since we made it this far.
TiXmlElement* mesh = docHandle.FirstChildElement("NEKTAR").FirstChildElement("GEOMETRY").Element();
TiXmlElement* field = NULL;
ASSERTL0(mesh, "Unable to find GEOMETRY tag in file.");
/// Look for elements in ELEMENT block.
field = mesh->FirstChildElement("COMPOSITE");
ASSERTL0(field, "Unable to find COMPOSITE tag in file.");
int nextCompositeNumber = -1;
/// All elements are of the form: "<C ID = "N"> ... </C>".
/// Read the ID field first.
TiXmlElement *composite = field->FirstChildElement("C");
while (composite)
{
nextCompositeNumber++;
int indx;
int err = composite->QueryIntAttribute("ID", &indx);
ASSERTL0(err == TIXML_SUCCESS, "Unable to read attribute ID.");
// ASSERTL0(indx == nextCompositeNumber, "Composite IDs must begin with zero and be sequential.");
TiXmlNode* compositeChild = composite->FirstChild();
// This is primarily to skip comments that may be present.
// Comments appear as nodes just like elements.
// We are specifically looking for text in the body
// of the definition.
while(compositeChild && compositeChild->Type() != TiXmlNode::TEXT)
{
compositeChild = compositeChild->NextSibling();
}
ASSERTL0(compositeChild, "Unable to read composite definition body.");
std::string compositeStr = compositeChild->ToText()->ValueStr();
/// Parse out the element components corresponding to type of element.
std::istringstream compositeDataStrm(compositeStr.c_str());
try
{
bool first = true;
std::string prevCompositeElementStr;
while (!compositeDataStrm.fail())
{
std::string compositeElementStr;
compositeDataStrm >> compositeElementStr;
if (!compositeDataStrm.fail())
{
if (first)
{
first = false;
Composite curVector(MemoryManager<GeometryVector>::AllocateSharedPtr());
m_meshComposites[indx] = curVector;
}
if (compositeElementStr.length() > 0)
{
ResolveGeomRef(prevCompositeElementStr, compositeElementStr, m_meshComposites[indx]);
}
prevCompositeElementStr = compositeElementStr;
}
}
}
catch(...)
{
(std::string("Unable to read COMPOSITE data for composite: ") + compositeStr).c_str());
}
/// Keep looking
composite = composite->NextSiblingElement("C");
}
}
void Nektar::SpatialDomains::MeshGraph3D::ReadEdges ( TiXmlDocument &  doc)
protected

We know we have it since we made it this far.

Look for elements in ELEMENT block.

All elements are of the form: "<E ID="#"> ... </E>", with ? being the element type. Read the ID field first.

Since all edge data is one big text block, we need to accumulate all TEXT data and then parse it. This approach effectively skips all comments or other node types since we only care about the edge list. We cannot handle missing edge numbers as we could with missing element numbers due to the text block format.

Now parse out the edges, three fields at a time.

Definition at line 98 of file MeshGraph3D.cpp.

References ASSERTL0, ErrorUtil::efatal, Nektar::SpatialDomains::MeshGraph::GetVertex(), Nektar::SpatialDomains::MeshGraph::m_curvedEdges, Nektar::SpatialDomains::MeshGraph::m_segGeoms, Nektar::SpatialDomains::MeshGraph::m_spaceDimension, and NEKERROR.

Referenced by ReadGeometry().

{
/// We know we have it since we made it this far.
TiXmlHandle docHandle(&doc);
TiXmlElement* mesh = docHandle.FirstChildElement("NEKTAR").FirstChildElement("GEOMETRY").Element();
TiXmlElement* field = NULL;
/// Look for elements in ELEMENT block.
field = mesh->FirstChildElement("EDGE");
ASSERTL0(field, "Unable to find EDGE tag in file.");
/// All elements are of the form: "<E ID="#"> ... </E>", with
/// ? being the element type.
/// Read the ID field first.
TiXmlElement *edge = field->FirstChildElement("E");
/// Since all edge data is one big text block, we need to accumulate
/// all TEXT data and then parse it. This approach effectively skips
/// all comments or other node types since we only care about the
/// edge list. We cannot handle missing edge numbers as we could
/// with missing element numbers due to the text block format.
std::string edgeStr;
int i,indx;
int nextEdgeNumber = -1;
// Curved Edges
map<int, int> edge_curved;
for(i = 0; i < m_curvedEdges.size(); ++i)
{
edge_curved[m_curvedEdges[i]->m_curveID] = i;
}
while(edge)
{
nextEdgeNumber++;
int err = edge->QueryIntAttribute("ID",&indx);
ASSERTL0(err == TIXML_SUCCESS, "Unable to read edge attribute ID.");
TiXmlNode *child = edge->FirstChild();
edgeStr.clear();
if (child->Type() == TiXmlNode::TEXT)
{
edgeStr += child->ToText()->ValueStr();
}
/// Now parse out the edges, three fields at a time.
int vertex1, vertex2;
std::istringstream edgeDataStrm(edgeStr.c_str());
try
{
while (!edgeDataStrm.fail())
{
edgeDataStrm >> vertex1 >> vertex2;
// Must check after the read because we may be at the end and not know it.
// If we are at the end we will add a duplicate of the last entry if we
// don't check here.
if (!edgeDataStrm.fail())
{
PointGeomSharedPtr vertices[2] = {GetVertex(vertex1), GetVertex(vertex2)};
if (edge_curved.count(indx) == 0)
{
edge = MemoryManager<SegGeom>::AllocateSharedPtr(indx, m_spaceDimension, vertices);
}
else
{
edge = MemoryManager<SegGeom>::AllocateSharedPtr(indx, m_spaceDimension, vertices, m_curvedEdges[edge_curved.find(indx)->second]);
}
m_segGeoms[indx] = edge;
}
}
}
catch(...)
{
NEKERROR(ErrorUtil::efatal, (std::string("Unable to read edge data: ") + edgeStr).c_str());
}
edge = edge->NextSiblingElement("E");
}
}
void Nektar::SpatialDomains::MeshGraph3D::ReadElements ( TiXmlDocument &  doc)
protected

We know we have it since we made it this far.

Look for elements in ELEMENT block.

All elements are of the form: "<? ID="#"> ... </?>", with ? being the element type.

These should be ordered.

Read id attribute.

Read text element description.

Parse out the element components corresponding to type of element.

Create arrays for the tri and quad faces.

Fill the arrays and make sure there aren't too many faces.

Make sure all of the face indicies could be read, and that there weren't too few.

Create arrays for the tri and quad faces.

Fill the arrays and make sure there aren't too many faces.

Make sure all of the face indicies could be read, and that there weren't too few.

Create arrays for the tri and quad faces.

Fill the arrays and make sure there aren't too many faces.

Make sure all of the face indicies could be read, and that there weren't too few.

Create arrays for the tri and quad faces.

Fill the arrays and make sure there aren't too many faces.

Make sure all of the face indicies could be read, and that there weren't too few.

Keep looking

Definition at line 338 of file MeshGraph3D.cpp.

References ASSERTL0, ErrorUtil::efatal, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eTriangle, GetGeometry2D(), Nektar::SpatialDomains::TetGeom::kNfaces, Nektar::SpatialDomains::PyrGeom::kNfaces, Nektar::SpatialDomains::PrismGeom::kNfaces, Nektar::SpatialDomains::HexGeom::kNfaces, Nektar::SpatialDomains::TetGeom::kNqfaces, Nektar::SpatialDomains::PyrGeom::kNqfaces, Nektar::SpatialDomains::PrismGeom::kNqfaces, Nektar::SpatialDomains::HexGeom::kNqfaces, Nektar::SpatialDomains::TetGeom::kNtfaces, Nektar::SpatialDomains::PyrGeom::kNtfaces, Nektar::SpatialDomains::PrismGeom::kNtfaces, Nektar::SpatialDomains::HexGeom::kNtfaces, Nektar::SpatialDomains::MeshGraph::m_hexGeoms, Nektar::SpatialDomains::MeshGraph::m_prismGeoms, Nektar::SpatialDomains::MeshGraph::m_pyrGeoms, Nektar::SpatialDomains::MeshGraph::m_tetGeoms, NEKERROR, and PopulateFaceToElMap().

Referenced by ReadGeometry().

{
/// We know we have it since we made it this far.
TiXmlHandle docHandle(&doc);
TiXmlElement* mesh = docHandle.FirstChildElement("NEKTAR").FirstChildElement("GEOMETRY").Element();
TiXmlElement* field = NULL;
/// Look for elements in ELEMENT block.
field = mesh->FirstChildElement("ELEMENT");
ASSERTL0(field, "Unable to find ELEMENT tag in file.");
int nextElementNumber = -1;
/// All elements are of the form: "<? ID="#"> ... </?>", with
/// ? being the element type.
TiXmlElement *element = field->FirstChildElement();
while (element)
{
std::string elementType(element->ValueStr());
//A - tet, P - pyramid, R - prism, H - hex
ASSERTL0(elementType == "A" || elementType == "P" || elementType == "R" || elementType == "H",
(std::string("Unknown 3D element type: ") + elementType).c_str());
/// These should be ordered.
nextElementNumber++;
/// Read id attribute.
int indx;
int err = element->QueryIntAttribute("ID", &indx);
ASSERTL0(err == TIXML_SUCCESS, "Unable to read element attribute ID.");
// ASSERTL0(indx == nextElementNumber, "Element IDs must begin with zero and be sequential.");
/// Read text element description.
TiXmlNode* elementChild = element->FirstChild();
std::string elementStr;
while(elementChild)
{
if (elementChild->Type() == TiXmlNode::TEXT)
{
elementStr += elementChild->ToText()->ValueStr();
}
elementChild = elementChild->NextSibling();
}
ASSERTL0(!elementStr.empty(), "Unable to read element description body.");
std::istringstream elementDataStrm(elementStr.c_str());
/// Parse out the element components corresponding to type of element.
// Tetrahedral
if (elementType == "A")
{
try
{
/// Create arrays for the tri and quad faces.
const int kNfaces = TetGeom::kNfaces;
const int kNtfaces = TetGeom::kNtfaces;
const int kNqfaces = TetGeom::kNqfaces;
TriGeomSharedPtr tfaces[kNtfaces];
//QuadGeomSharedPtr qfaces[kNqfaces];
int Ntfaces = 0;
int Nqfaces = 0;
/// Fill the arrays and make sure there aren't too many faces.
std::stringstream errorstring;
errorstring << "Element " << indx << " must have " << kNtfaces << " triangle face(s), and " << kNqfaces << " quadrilateral face(s).";
for (int i = 0; i < kNfaces; i++)
{
int faceID;
elementDataStrm >> faceID;
if (face == Geometry2DSharedPtr() ||
(face->GetShapeType() != LibUtilities::eTriangle && face->GetShapeType() != LibUtilities::eQuadrilateral))
{
std::stringstream errorstring;
errorstring << "Element " << indx << " has invalid face: " << faceID;
ASSERTL0(false, errorstring.str().c_str());
}
else if (face->GetShapeType() == LibUtilities::eTriangle)
{
ASSERTL0(Ntfaces < kNtfaces, errorstring.str().c_str());
tfaces[Ntfaces++] = boost::static_pointer_cast<TriGeom>(face);
}
else if (face->GetShapeType() == LibUtilities::eQuadrilateral)
{
ASSERTL0(Nqfaces < kNqfaces, errorstring.str().c_str());
}
}
/// Make sure all of the face indicies could be read, and that there weren't too few.
ASSERTL0(!elementDataStrm.fail(), (std::string("Unable to read element data for TETRAHEDRON: ") + elementStr).c_str());
ASSERTL0(Ntfaces == kNtfaces, errorstring.str().c_str());
ASSERTL0(Nqfaces == kNqfaces, errorstring.str().c_str());
TetGeomSharedPtr tetgeom(MemoryManager<TetGeom>::AllocateSharedPtr(tfaces));
tetgeom->SetGlobalID(indx);
m_tetGeoms[indx] = tetgeom;
PopulateFaceToElMap(tetgeom, kNfaces);
}
catch(...)
{
(std::string("Unable to read element data for TETRAHEDRON: ") + elementStr).c_str());
}
}
// Pyramid
else if (elementType == "P")
{
try
{
/// Create arrays for the tri and quad faces.
const int kNfaces = PyrGeom::kNfaces;
const int kNtfaces = PyrGeom::kNtfaces;
const int kNqfaces = PyrGeom::kNqfaces;
Geometry2DSharedPtr faces[kNfaces];
int Nfaces = 0;
int Ntfaces = 0;
int Nqfaces = 0;
/// Fill the arrays and make sure there aren't too many faces.
std::stringstream errorstring;
errorstring << "Element " << indx << " must have " << kNtfaces << " triangle face(s), and " << kNqfaces << " quadrilateral face(s).";
for (int i = 0; i < kNfaces; i++)
{
int faceID;
elementDataStrm >> faceID;
if (face == Geometry2DSharedPtr() ||
(face->GetShapeType() != LibUtilities::eTriangle && face->GetShapeType() != LibUtilities::eQuadrilateral))
{
std::stringstream errorstring;
errorstring << "Element " << indx << " has invalid face: " << faceID;
ASSERTL0(false, errorstring.str().c_str());
}
else if (face->GetShapeType() == LibUtilities::eTriangle)
{
ASSERTL0(Ntfaces < kNtfaces, errorstring.str().c_str());
faces[Nfaces++] = boost::static_pointer_cast<TriGeom>(face);
Ntfaces++;
}
else if (face->GetShapeType() == LibUtilities::eQuadrilateral)
{
ASSERTL0(Nqfaces < kNqfaces, errorstring.str().c_str());
faces[Nfaces++] = boost::static_pointer_cast<QuadGeom>(face);
Nqfaces++;
}
}
/// Make sure all of the face indicies could be read, and that there weren't too few.
ASSERTL0(!elementDataStrm.fail(), (std::string("Unable to read element data for PYRAMID: ") + elementStr).c_str());
ASSERTL0(Ntfaces == kNtfaces, errorstring.str().c_str());
ASSERTL0(Nqfaces == kNqfaces, errorstring.str().c_str());
PyrGeomSharedPtr pyrgeom(MemoryManager<PyrGeom>::AllocateSharedPtr(faces));
pyrgeom->SetGlobalID(indx);
m_pyrGeoms[indx] = pyrgeom;
PopulateFaceToElMap(pyrgeom, kNfaces);
}
catch(...)
{
(std::string("Unable to read element data for PYRAMID: ") + elementStr).c_str());
}
}
// Prism
else if (elementType == "R")
{
try
{
/// Create arrays for the tri and quad faces.
const int kNfaces = PrismGeom::kNfaces;
const int kNtfaces = PrismGeom::kNtfaces;
const int kNqfaces = PrismGeom::kNqfaces;
Geometry2DSharedPtr faces[kNfaces];
int Ntfaces = 0;
int Nqfaces = 0;
int Nfaces = 0;
/// Fill the arrays and make sure there aren't too many faces.
std::stringstream errorstring;
errorstring << "Element " << indx << " must have "
<< kNtfaces << " triangle face(s), and "
<< kNqfaces << " quadrilateral face(s).";
for (int i = 0; i < kNfaces; i++)
{
int faceID;
elementDataStrm >> faceID;
if (face == Geometry2DSharedPtr() ||
(face->GetShapeType() != LibUtilities::eTriangle && face->GetShapeType() != LibUtilities::eQuadrilateral))
{
std::stringstream errorstring;
errorstring << "Element " << indx << " has invalid face: " << faceID;
ASSERTL0(false, errorstring.str().c_str());
}
else if (face->GetShapeType() == LibUtilities::eTriangle)
{
ASSERTL0(Ntfaces < kNtfaces, errorstring.str().c_str());
faces[Nfaces++] = boost::static_pointer_cast<TriGeom>(face);
Ntfaces++;
}
else if (face->GetShapeType() == LibUtilities::eQuadrilateral)
{
ASSERTL0(Nqfaces < kNqfaces, errorstring.str().c_str());
faces[Nfaces++] = boost::static_pointer_cast<QuadGeom>(face);
Nqfaces++;
}
}
/// Make sure all of the face indicies could be read, and that there weren't too few.
ASSERTL0(!elementDataStrm.fail(), (std::string("Unable to read element data for PRISM: ") + elementStr).c_str());
ASSERTL0(Ntfaces == kNtfaces, errorstring.str().c_str());
ASSERTL0(Nqfaces == kNqfaces, errorstring.str().c_str());
PrismGeomSharedPtr prismgeom(MemoryManager<PrismGeom>::AllocateSharedPtr(faces));
prismgeom->SetGlobalID(indx);
m_prismGeoms[indx] = prismgeom;
PopulateFaceToElMap(prismgeom, kNfaces);
}
catch(...)
{
(std::string("Unable to read element data for PRISM: ") + elementStr).c_str());
}
}
// Hexahedral
else if (elementType == "H")
{
try
{
/// Create arrays for the tri and quad faces.
const int kNfaces = HexGeom::kNfaces;
const int kNtfaces = HexGeom::kNtfaces;
const int kNqfaces = HexGeom::kNqfaces;
//TriGeomSharedPtr tfaces[kNtfaces];
QuadGeomSharedPtr qfaces[kNqfaces];
int Ntfaces = 0;
int Nqfaces = 0;
/// Fill the arrays and make sure there aren't too many faces.
std::stringstream errorstring;
errorstring << "Element " << indx << " must have " << kNtfaces << " triangle face(s), and " << kNqfaces << " quadrilateral face(s).";
for (int i = 0; i < kNfaces; i++)
{
int faceID;
elementDataStrm >> faceID;
if (face == Geometry2DSharedPtr() ||
(face->GetShapeType() != LibUtilities::eTriangle && face->GetShapeType() != LibUtilities::eQuadrilateral))
{
std::stringstream errorstring;
errorstring << "Element " << indx << " has invalid face: " << faceID;
ASSERTL0(false, errorstring.str().c_str());
}
else if (face->GetShapeType() == LibUtilities::eTriangle)
{
ASSERTL0(Ntfaces < kNtfaces, errorstring.str().c_str());
//tfaces[Ntfaces++] = boost::static_pointer_cast<TriGeom>(face);
}
else if (face->GetShapeType() == LibUtilities::eQuadrilateral)
{
ASSERTL0(Nqfaces < kNqfaces, errorstring.str().c_str());
qfaces[Nqfaces++] = boost::static_pointer_cast<QuadGeom>(face);
}
}
/// Make sure all of the face indicies could be read, and that there weren't too few.
ASSERTL0(!elementDataStrm.fail(), (std::string("Unable to read element data for HEXAHEDRAL: ") + elementStr).c_str());
ASSERTL0(Ntfaces == kNtfaces, errorstring.str().c_str());
ASSERTL0(Nqfaces == kNqfaces, errorstring.str().c_str());
HexGeomSharedPtr hexgeom(MemoryManager<HexGeom>::AllocateSharedPtr(qfaces));
hexgeom->SetGlobalID(indx);
m_hexGeoms[indx] = hexgeom;
PopulateFaceToElMap(hexgeom, kNfaces);
}
catch(...)
{
(std::string("Unable to read element data for HEXAHEDRAL: ") + elementStr).c_str());
}
}
/// Keep looking
element = element->NextSiblingElement();
}
}
void Nektar::SpatialDomains::MeshGraph3D::ReadFaces ( TiXmlDocument &  doc)
protected

We know we have it since we made it this far.

Look for elements in FACE block.

All faces are of the form: "<? ID="#"> ... </?>", with ? being an element type (either Q or T).

Read id attribute.

Read text element description.

Parse out the element components corresponding to type of element.

Create a TriGeom to hold the new definition.

Create a QuadGeom to hold the new definition.

Keep looking

Definition at line 185 of file MeshGraph3D.cpp.

References ASSERTL0, ErrorUtil::efatal, Nektar::SpatialDomains::SegGeom::GetEdgeOrientation(), GetSegGeom(), Nektar::SpatialDomains::TriGeom::kNedges, Nektar::SpatialDomains::QuadGeom::kNedges, Nektar::SpatialDomains::MeshGraph::m_curvedFaces, Nektar::SpatialDomains::MeshGraph::m_quadGeoms, Nektar::SpatialDomains::MeshGraph::m_triGeoms, and NEKERROR.

Referenced by ReadGeometry().

{
/// We know we have it since we made it this far.
TiXmlHandle docHandle(&doc);
TiXmlElement* mesh = docHandle.FirstChildElement("NEKTAR").FirstChildElement("GEOMETRY").Element();
TiXmlElement* field = NULL;
/// Look for elements in FACE block.
field = mesh->FirstChildElement("FACE");
ASSERTL0(field, "Unable to find FACE tag in file.");
// Curved faces
map<int, int> face_curved;
for (int i = 0; i < m_curvedFaces.size(); ++i)
{
face_curved[m_curvedFaces[i]->m_curveID] = i;
}
/// All faces are of the form: "<? ID="#"> ... </?>", with
/// ? being an element type (either Q or T).
TiXmlElement *element = field->FirstChildElement();
while (element)
{
std::string elementType(element->ValueStr());
ASSERTL0(elementType == "Q" || elementType == "T",
(std::string("Unknown 3D face type: ") + elementType).c_str());
/// Read id attribute.
int indx;
int err = element->QueryIntAttribute("ID", &indx);
ASSERTL0(err == TIXML_SUCCESS, "Unable to read face attribute ID.");
/// Read text element description.
TiXmlNode* elementChild = element->FirstChild();
std::string elementStr;
while(elementChild)
{
if (elementChild->Type() == TiXmlNode::TEXT)
{
elementStr += elementChild->ToText()->ValueStr();
}
elementChild = elementChild->NextSibling();
}
ASSERTL0(!elementStr.empty(), "Unable to read face description body.");
/// Parse out the element components corresponding to type of element.
if (elementType == "T")
{
// Read three edge numbers
int edge1, edge2, edge3;
std::istringstream elementDataStrm(elementStr.c_str());
try
{
elementDataStrm >> edge1;
elementDataStrm >> edge2;
elementDataStrm >> edge3;
ASSERTL0(!elementDataStrm.fail(), (std::string("Unable to read face data for TRIANGLE: ") + elementStr).c_str());
/// Create a TriGeom to hold the new definition.
{
GetSegGeom(edge1),
GetSegGeom(edge2),
GetSegGeom(edge3)
};
{
SegGeom::GetEdgeOrientation(*edges[0], *edges[1]),
SegGeom::GetEdgeOrientation(*edges[1], *edges[2]),
SegGeom::GetEdgeOrientation(*edges[2], *edges[0])
};
if (face_curved.count(indx) == 0)
{
trigeom = MemoryManager<TriGeom>::AllocateSharedPtr(indx, edges, edgeorient);
}
else
{
trigeom = MemoryManager<TriGeom>::AllocateSharedPtr(indx, edges, edgeorient, m_curvedFaces[face_curved.find(indx)->second]);
}
trigeom->SetGlobalID(indx);
m_triGeoms[indx] = trigeom;
}
catch(...)
{
(std::string("Unable to read face data for TRIANGLE: ") + elementStr).c_str());
}
}
else if (elementType == "Q")
{
// Read four edge numbers
int edge1, edge2, edge3, edge4;
std::istringstream elementDataStrm(elementStr.c_str());
try
{
elementDataStrm >> edge1;
elementDataStrm >> edge2;
elementDataStrm >> edge3;
elementDataStrm >> edge4;
ASSERTL0(!elementDataStrm.fail(), (std::string("Unable to read face data for QUAD: ") + elementStr).c_str());
/// Create a QuadGeom to hold the new definition.
{GetSegGeom(edge1),GetSegGeom(edge2),
GetSegGeom(edge3),GetSegGeom(edge4)};
{
SegGeom::GetEdgeOrientation(*edges[0], *edges[1]),
SegGeom::GetEdgeOrientation(*edges[1], *edges[2]),
SegGeom::GetEdgeOrientation(*edges[2], *edges[3]),
SegGeom::GetEdgeOrientation(*edges[3], *edges[0])
};
if (face_curved.count(indx) == 0)
{
quadgeom = MemoryManager<QuadGeom>::AllocateSharedPtr(indx, edges, edgeorient);
} else {
quadgeom = MemoryManager<QuadGeom>::AllocateSharedPtr(indx, edges, edgeorient, m_curvedFaces[face_curved.find(indx)->second]);
}
quadgeom->SetGlobalID(indx);
m_quadGeoms[indx] = quadgeom;
}
catch(...)
{
NEKERROR(ErrorUtil::efatal,(std::string("Unable to read face data for QUAD: ") + elementStr).c_str());
}
}
/// Keep looking
element = element->NextSiblingElement();
}
}
void Nektar::SpatialDomains::MeshGraph3D::ReadGeometry ( const std::string &  infilename)
virtual

Read will read the meshgraph vertices given a filename.

Reimplemented from Nektar::SpatialDomains::MeshGraph.

Definition at line 62 of file MeshGraph3D.cpp.

References ASSERTL0.

Referenced by MeshGraph3D(), and ReadGeometry().

{
TiXmlDocument doc(infilename);
bool loadOkay = doc.LoadFile();
std::stringstream errstr;
errstr << "Unable to load file: " << infilename << "\n";
errstr << doc.ErrorDesc() << " (Line " << doc.ErrorRow()
<< ", Column " << doc.ErrorCol() << ")";
ASSERTL0(loadOkay, errstr.str());
}
void Nektar::SpatialDomains::MeshGraph3D::ReadGeometry ( TiXmlDocument &  doc)
virtual

Read will read the meshgraph vertices given a TiXmlDocument.

Look for all geometry related data in GEOMETRY block.

Error value returned by TinyXML.

Reimplemented from Nektar::SpatialDomains::MeshGraph.

Definition at line 77 of file MeshGraph3D.cpp.

References ASSERTL0, ReadComposites(), Nektar::SpatialDomains::MeshGraph::ReadCurves(), Nektar::SpatialDomains::MeshGraph::ReadDomain(), ReadEdges(), ReadElements(), ReadFaces(), and ReadGeometry().

{
// Read mesh first
TiXmlHandle docHandle(&doc);
TiXmlElement* mesh = NULL;
/// Look for all geometry related data in GEOMETRY block.
mesh = docHandle.FirstChildElement("NEKTAR").FirstChildElement("GEOMETRY").Element();
ASSERTL0(mesh, "Unable to find GEOMETRY tag in file.");
ReadCurves(doc);
ReadEdges(doc);
ReadFaces(doc);
ReadDomain(doc);
}
void Nektar::SpatialDomains::MeshGraph3D::ResolveGeomRef ( const std::string &  prevToken,
const std::string &  token,
Composite composite 
)
protected

Definition at line 755 of file MeshGraph3D.cpp.

References ASSERTL0, Nektar::SpatialDomains::MeshGraph::CheckRange(), ErrorUtil::efatal, ErrorUtil::ewarning, Nektar::ParseUtils::GenerateSeqVector(), GetGeometry2D(), Nektar::iterator, Nektar::SpatialDomains::MeshGraph::m_hexGeoms, Nektar::SpatialDomains::MeshGraph::m_prismGeoms, Nektar::SpatialDomains::MeshGraph::m_pyrGeoms, Nektar::SpatialDomains::MeshGraph::m_quadGeoms, Nektar::SpatialDomains::MeshGraph::m_segGeoms, Nektar::SpatialDomains::MeshGraph::m_tetGeoms, Nektar::SpatialDomains::MeshGraph::m_triGeoms, Nektar::SpatialDomains::MeshGraph::m_vertSet, and NEKERROR.

Referenced by ReadComposites().

{
try
{
std::istringstream tokenStream(token);
std::istringstream prevTokenStream(prevToken);
char type;
char prevType;
tokenStream >> type;
std::string::size_type indxBeg = token.find_first_of('[') + 1;
std::string::size_type indxEnd = token.find_last_of(']') - 1;
ASSERTL0(indxBeg <= indxEnd, (std::string("Error reading index definition:") + token).c_str());
std::string indxStr = token.substr(indxBeg, indxEnd - indxBeg + 1);
std::vector<unsigned int> seqVector;
bool err = ParseUtils::GenerateSeqVector(indxStr.c_str(), seqVector);
ASSERTL0(err, (std::string("Error reading composite elements: ") + indxStr).c_str());
prevTokenStream >> prevType;
// All composites must be of the same dimension. This map makes things clean to compare.
map<char, int> typeMap;
typeMap['V'] = 1; // Vertex
typeMap['E'] = 1; // Edge
typeMap['T'] = 2; // Triangle
typeMap['Q'] = 2; // Quad
typeMap['A'] = 3; // Tet
typeMap['P'] = 3; // Pyramid
typeMap['R'] = 3; // Prism
typeMap['H'] = 3; // Hex
// Make sure only geoms of the same dimension are combined.
bool validSequence = (prevToken.empty() || (typeMap[type] == typeMap[prevType]));
ASSERTL0(validSequence, (std::string("Invalid combination of composite items: ")
+ type + " and " + prevType + ".").c_str());
switch(type)
{
case 'V': // Vertex
for (seqIter = seqVector.begin(); seqIter != seqVector.end(); ++seqIter)
{
if (m_vertSet.find(*seqIter) == m_vertSet.end())
{
char errStr[16] = "";
::sprintf(errStr, "%d", *seqIter);
NEKERROR(ErrorUtil::ewarning, (std::string("Unknown vertex index: ") + errStr).c_str());
}
else
{
composite->push_back(m_vertSet[*seqIter]);
}
}
break;
case 'E': // Edge
for (seqIter = seqVector.begin(); seqIter != seqVector.end(); ++seqIter)
{
if (m_segGeoms.find(*seqIter) == m_segGeoms.end())
{
char errStr[16] = "";
::sprintf(errStr, "%d", *seqIter);
NEKERROR(ErrorUtil::ewarning, (std::string("Unknown edge index: ") + errStr).c_str());
}
else
{
composite->push_back(m_segGeoms[*seqIter]);
}
}
break;
case 'F': // Face
for (seqIter = seqVector.begin(); seqIter != seqVector.end(); ++seqIter)
{
if (face == Geometry2DSharedPtr())
{
char errStr[16] = "";
::sprintf(errStr, "%d", *seqIter);
NEKERROR(ErrorUtil::ewarning, (std::string("Unknown face index: ") + errStr).c_str());
}
else
{
if(CheckRange(*face))
{
composite->push_back(face);
}
}
}
break;
case 'T': // Triangle
for (seqIter = seqVector.begin(); seqIter != seqVector.end(); ++seqIter)
{
if (m_triGeoms.find(*seqIter) == m_triGeoms.end())
{
char errStr[16] = "";
::sprintf(errStr, "%d", *seqIter);
NEKERROR(ErrorUtil::ewarning, (std::string("Unknown triangle index: ") + errStr).c_str());
}
else
{
if(CheckRange(*m_triGeoms[*seqIter]))
{
composite->push_back(m_triGeoms[*seqIter]);
}
}
}
break;
case 'Q': // Quad
for (seqIter = seqVector.begin(); seqIter != seqVector.end(); ++seqIter)
{
if (m_quadGeoms.find(*seqIter) == m_quadGeoms.end())
{
char errStr[16] = "";
::sprintf(errStr, "%d", *seqIter);
NEKERROR(ErrorUtil::ewarning, (std::string("Unknown quad index: ") + errStr).c_str());
}
else
{
if(CheckRange(*m_quadGeoms[*seqIter]))
{
composite->push_back(m_quadGeoms[*seqIter]);
}
}
}
break;
// Tetrahedron
case 'A':
for (seqIter = seqVector.begin(); seqIter != seqVector.end(); ++seqIter)
{
if (m_tetGeoms.find(*seqIter) == m_tetGeoms.end())
{
char errStr[16] = "";
::sprintf(errStr, "%d", *seqIter);
NEKERROR(ErrorUtil::ewarning, (std::string("Unknown tet index: ") + errStr).c_str());
}
else
{
if(CheckRange(*m_tetGeoms[*seqIter]))
{
composite->push_back(m_tetGeoms[*seqIter]);
}
}
}
break;
// Pyramid
case 'P':
for (seqIter = seqVector.begin(); seqIter != seqVector.end(); ++seqIter)
{
if (m_pyrGeoms.find(*seqIter) == m_pyrGeoms.end())
{
char errStr[16] = "";
::sprintf(errStr, "%d", *seqIter);
NEKERROR(ErrorUtil::ewarning, (std::string("Unknown pyramid index: ") + errStr).c_str());
}
else
{
if(CheckRange(*m_pyrGeoms[*seqIter]))
{
composite->push_back(m_pyrGeoms[*seqIter]);
}
}
}
break;
// Prism
case 'R':
for (seqIter = seqVector.begin(); seqIter != seqVector.end(); ++seqIter)
{
if (m_prismGeoms.find(*seqIter) == m_prismGeoms.end())
{
char errStr[16] = "";
::sprintf(errStr, "%d", *seqIter);
NEKERROR(ErrorUtil::ewarning, (std::string("Unknown prism index: ") + errStr).c_str());
}
else
{
if(CheckRange(*m_prismGeoms[*seqIter]))
{
composite->push_back(m_prismGeoms[*seqIter]);
}
}
}
break;
// Hex
case 'H':
for (seqIter = seqVector.begin(); seqIter != seqVector.end(); ++seqIter)
{
if (m_hexGeoms.find(*seqIter) == m_hexGeoms.end())
{
char errStr[16] = "";
::sprintf(errStr, "%d", *seqIter);
NEKERROR(ErrorUtil::ewarning, (std::string("Unknown hex index: ") + errStr).c_str());
}
else
{
if(CheckRange(*m_hexGeoms[*seqIter]))
{
composite->push_back(m_hexGeoms[*seqIter]);
}
}
}
break;
default:
NEKERROR(ErrorUtil::efatal, (std::string("Unrecognized composite token: ") + token).c_str());
}
}
catch(...)
{
NEKERROR(ErrorUtil::efatal, (std::string("Problem processing composite token: ") + token).c_str());
}
return;
}

Member Data Documentation

boost::unordered_map<int, ElementFaceVectorSharedPtr> Nektar::SpatialDomains::MeshGraph3D::m_faceToElMap
private

Definition at line 220 of file MeshGraph3D.h.

Referenced by GetElementsFromFace(), and PopulateFaceToElMap().