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

Base class for a spectral/hp element mesh. More...

#include <MeshGraph.h>

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

Public Member Functions

 MeshGraph ()
 MeshGraph (unsigned int meshDimension, unsigned int spaceDimension)
 MeshGraph (const LibUtilities::SessionReaderSharedPtr &pSession, const DomainRangeShPtr &rng=NullDomainRangeShPtr)
virtual ~MeshGraph ()
virtual void ReadGeometry (const std::string &infilename)
 Read will read the meshgraph vertices given a filename.
virtual void ReadGeometry (TiXmlDocument &doc)
 Read will read the meshgraph vertices given a TiXmlDocument.
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.

Static Public Member Functions

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 Member Functions

ExpansionMapShPtr SetUpExpansionMap (void)

Protected Attributes

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

Base class for a spectral/hp element mesh.

Definition at line 180 of file MeshGraph.h.

Constructor & Destructor Documentation

Nektar::SpatialDomains::MeshGraph::MeshGraph ( )
Nektar::SpatialDomains::MeshGraph::MeshGraph ( unsigned int  meshDimension,
unsigned int  spaceDimension 
)

Definition at line 84 of file MeshGraph.cpp.

:
m_meshDimension(meshDimension),
m_spaceDimension(spaceDimension),
{
}
Nektar::SpatialDomains::MeshGraph::MeshGraph ( const LibUtilities::SessionReaderSharedPtr pSession,
const DomainRangeShPtr rng = NullDomainRangeShPtr 
)

Definition at line 97 of file MeshGraph.cpp.

:
m_session(pSession),
{
}
Nektar::SpatialDomains::MeshGraph::~MeshGraph ( )
virtual

Definition at line 110 of file MeshGraph.cpp.

{
}

Member Function Documentation

SegGeomSharedPtr Nektar::SpatialDomains::MeshGraph::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.

Definition at line 3386 of file MeshGraph.cpp.

References m_segGeoms, and m_spaceDimension.

{
PointGeomSharedPtr vertices[] = {v0, v1};
int edgeId = m_segGeoms.rbegin()->first + 1;
if( curveDefinition )
{
edge = MemoryManager<SegGeom>::AllocateSharedPtr(edgeId, m_spaceDimension, vertices, curveDefinition);
}
else
{
edge = MemoryManager<SegGeom>::AllocateSharedPtr(edgeId, m_spaceDimension, vertices);
}
m_segGeoms[edgeId] = edge;
return edge;
}
HexGeomSharedPtr Nektar::SpatialDomains::MeshGraph::AddHexahedron ( QuadGeomSharedPtr  qfaces[HexGeom::kNqfaces])

Definition at line 3486 of file MeshGraph.cpp.

References m_hexGeoms.

{
unsigned int index = m_hexGeoms.rbegin()->first + 1;
HexGeomSharedPtr hexgeom(MemoryManager<HexGeom>::AllocateSharedPtr(qfaces));
hexgeom->SetGlobalID(index);
m_hexGeoms[index] = hexgeom;
return hexgeom;
}
PrismGeomSharedPtr Nektar::SpatialDomains::MeshGraph::AddPrism ( TriGeomSharedPtr  tfaces[PrismGeom::kNtfaces],
QuadGeomSharedPtr  qfaces[PrismGeom::kNqfaces] 
)

Definition at line 3438 of file MeshGraph.cpp.

References m_prismGeoms.

{
// Setting the orientation is disabled in the reader. Why?
Geometry2DSharedPtr faces[] = { qfaces[0], tfaces[0], qfaces[1], tfaces[1], qfaces[2] };
unsigned int index = m_prismGeoms.rbegin()->first + 1;
PrismGeomSharedPtr prismgeom(MemoryManager<PrismGeom>::AllocateSharedPtr(faces));
prismgeom->SetGlobalID(index);
m_prismGeoms[index] = prismgeom;
return prismgeom;
}
PyrGeomSharedPtr Nektar::SpatialDomains::MeshGraph::AddPyramid ( TriGeomSharedPtr  tfaces[PyrGeom::kNtfaces],
QuadGeomSharedPtr  qfaces[PyrGeom::kNqfaces] 
)

Definition at line 3469 of file MeshGraph.cpp.

References m_pyrGeoms.

{
Geometry2DSharedPtr faces[] = { qfaces[0], tfaces[0], tfaces[1], tfaces[2], tfaces[3] };
unsigned int index = m_pyrGeoms.rbegin()->first + 1;
PyrGeomSharedPtr pyrgeom(MemoryManager<PyrGeom>::AllocateSharedPtr(faces));
pyrgeom->SetGlobalID(index);
m_pyrGeoms[index] = pyrgeom;
return pyrgeom;
}
QuadGeomSharedPtr Nektar::SpatialDomains::MeshGraph::AddQuadrilateral ( SegGeomSharedPtr  edges[],
StdRegions::Orientation  orient[] 
)

Definition at line 3424 of file MeshGraph.cpp.

References m_quadGeoms.

{
int indx = m_quadGeoms.rbegin()->first + 1;
QuadGeomSharedPtr quadgeom(MemoryManager<QuadGeom>::AllocateSharedPtr(indx, edges, orient));
quadgeom->SetGlobalID(indx);
m_quadGeoms[indx] = quadgeom;
return quadgeom;
}
TetGeomSharedPtr Nektar::SpatialDomains::MeshGraph::AddTetrahedron ( TriGeomSharedPtr  tfaces[TetGeom::kNtfaces])

Definition at line 3455 of file MeshGraph.cpp.

References m_tetGeoms.

{
unsigned int index = m_tetGeoms.rbegin()->first + 1;
TetGeomSharedPtr tetgeom(MemoryManager<TetGeom>::AllocateSharedPtr(tfaces));
tetgeom->SetGlobalID(index);
m_tetGeoms[index] = tetgeom;
return tetgeom;
}
TriGeomSharedPtr Nektar::SpatialDomains::MeshGraph::AddTriangle ( SegGeomSharedPtr  edges[],
StdRegions::Orientation  orient[] 
)

Definition at line 3409 of file MeshGraph.cpp.

References m_triGeoms.

{
int indx = m_triGeoms.rbegin()->first + 1;
TriGeomSharedPtr trigeom(MemoryManager<TriGeom>::AllocateSharedPtr(indx, edges, orient));
trigeom->SetGlobalID(indx);
m_triGeoms[indx] = trigeom;
return trigeom;
}
PointGeomSharedPtr Nektar::SpatialDomains::MeshGraph::AddVertex ( NekDouble  x,
NekDouble  y,
NekDouble  z 
)

Adds a vertex to the with the next available ID.

Definition at line 3374 of file MeshGraph.cpp.

References m_spaceDimension, and m_vertSet.

{
unsigned int nextId = m_vertSet.rbegin()->first + 1;
PointGeomSharedPtr vert(MemoryManager<PointGeom>::AllocateSharedPtr(m_spaceDimension, nextId, x, y, z));
m_vertSet[nextId] = vert;
return vert;
}
bool Nektar::SpatialDomains::MeshGraph::CheckForGeomInfo ( std::string  parameter)
inline

Definition at line 523 of file MeshGraph.h.

References m_geomInfo.

{
return m_geomInfo.find(parameter) != m_geomInfo.end();
}
bool Nektar::SpatialDomains::MeshGraph::CheckRange ( Geometry2D geom)

Check if goemetry is in range definition if activated.

Definition at line 1415 of file MeshGraph.cpp.

References Nektar::SpatialDomains::Geometry::GetCoordim(), Nektar::SpatialDomains::Geometry::GetNumVerts(), Nektar::SpatialDomains::Geometry::GetVertex(), m_domainRange, and Nektar::SpatialDomains::NullDomainRangeShPtr.

Referenced by Nektar::SpatialDomains::MeshGraph3D::ResolveGeomRef(), and Nektar::SpatialDomains::MeshGraph2D::ResolveGeomRef().

{
bool returnval = true;
{
int nverts = geom.GetNumVerts();
int coordim = geom.GetCoordim();
// exclude elements outside x range if all vertices not in region
if(m_domainRange->doXrange)
{
int ncnt_low = 0;
int ncnt_up = 0;
for(int i = 0; i < nverts; ++i)
{
NekDouble xval = (*geom.GetVertex(i))[0];
if(xval < m_domainRange->xmin)
{
ncnt_low++;
}
if(xval > m_domainRange->xmax)
{
ncnt_up++;
}
}
// check for all verts to be less or greater than
// range so that if element spans thin range then
// it is still included
if((ncnt_up == nverts)||(ncnt_low == nverts))
{
returnval = false;
}
}
// exclude elements outside y range if all vertices not in region
if(m_domainRange->doYrange)
{
int ncnt_low = 0;
int ncnt_up = 0;
for(int i = 0; i < nverts; ++i)
{
NekDouble yval = (*geom.GetVertex(i))[1];
if(yval < m_domainRange->ymin)
{
ncnt_low++;
}
if(yval > m_domainRange->ymax)
{
ncnt_up++;
}
}
// check for all verts to be less or greater than
// range so that if element spans thin range then
// it is still included
if((ncnt_up == nverts)||(ncnt_low == nverts))
{
returnval = false;
}
}
if(coordim > 2)
{
// exclude elements outside z range if all vertices not in region
if(m_domainRange->doZrange)
{
int ncnt_low = 0;
int ncnt_up = 0;
for(int i = 0; i < nverts; ++i)
{
NekDouble zval = (*geom.GetVertex(i))[2];
if(zval < m_domainRange->zmin)
{
ncnt_low++;
}
if(zval > m_domainRange->zmax)
{
ncnt_up++;
}
}
// check for all verts to be less or greater than
// range so that if element spans thin range then
// it is still included
if((ncnt_up == nverts)||(ncnt_low == nverts))
{
returnval = false;
}
}
}
}
return returnval;
}
bool Nektar::SpatialDomains::MeshGraph::CheckRange ( Geometry3D geom)

Check if goemetry is in range definition if activated.

Definition at line 1518 of file MeshGraph.cpp.

References Nektar::SpatialDomains::Geometry::GetNumVerts(), Nektar::SpatialDomains::Geometry::GetVertex(), m_domainRange, and Nektar::SpatialDomains::NullDomainRangeShPtr.

{
bool returnval = true;
{
int nverts = geom.GetNumVerts();
if(m_domainRange->doXrange)
{
int ncnt_low = 0;
int ncnt_up = 0;
for(int i = 0; i < nverts; ++i)
{
NekDouble xval = (*geom.GetVertex(i))[0];
if(xval < m_domainRange->xmin)
{
ncnt_low++;
}
if(xval > m_domainRange->xmax)
{
ncnt_up++;
}
}
// check for all verts to be less or greater than
// range so that if element spans thin range then
// it is still included
if((ncnt_up == nverts)||(ncnt_low == nverts))
{
returnval = false;
}
}
if(m_domainRange->doYrange)
{
int ncnt_low = 0;
int ncnt_up = 0;
for(int i = 0; i < nverts; ++i)
{
NekDouble yval = (*geom.GetVertex(i))[1];
if(yval < m_domainRange->ymin)
{
ncnt_low++;
}
if(yval > m_domainRange->ymax)
{
ncnt_up++;
}
}
// check for all verts to be less or greater than
// range so that if element spans thin range then
// it is still included
if((ncnt_up == nverts)||(ncnt_low == nverts))
{
returnval = false;
}
}
if(m_domainRange->doZrange)
{
int ncnt_low = 0;
int ncnt_up = 0;
for(int i = 0; i < nverts; ++i)
{
NekDouble zval = (*geom.GetVertex(i))[2];
if(zval < m_domainRange->zmin)
{
ncnt_low++;
}
if(zval > m_domainRange->zmax)
{
ncnt_up++;
}
}
// check for all verts to be less or greater than
// range so that if element spans thin range then
// it is still included
if((ncnt_up == nverts)||(ncnt_low == nverts))
{
returnval = false;
}
}
}
return returnval;
}
LibUtilities::BasisKeyVector Nektar::SpatialDomains::MeshGraph::DefineBasisKeyFromExpansionType ( GeometrySharedPtr  in,
ExpansionType  type,
const int  order 
)
static

Definition at line 2586 of file MeshGraph.cpp.

References ASSERTL0, Nektar::LibUtilities::eChebyshev, Nektar::SpatialDomains::eChebyshev, Nektar::SpatialDomains::eChebyshevFourier, Nektar::LibUtilities::eFourier, Nektar::SpatialDomains::eFourier, Nektar::SpatialDomains::eFourierChebyshev, Nektar::LibUtilities::eFourierEvenlySpaced, Nektar::LibUtilities::eFourierHalfModeIm, Nektar::SpatialDomains::eFourierHalfModeIm, Nektar::LibUtilities::eFourierHalfModeRe, Nektar::SpatialDomains::eFourierHalfModeRe, Nektar::SpatialDomains::eFourierModified, Nektar::LibUtilities::eFourierSingleMode, Nektar::SpatialDomains::eFourierSingleMode, Nektar::LibUtilities::eFourierSingleModeSpaced, Nektar::LibUtilities::eGauss_Lagrange, Nektar::SpatialDomains::eGauss_Lagrange, Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eGaussGaussLegendre, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauMAlpha2Beta0, Nektar::LibUtilities::eGLL_Lagrange, Nektar::SpatialDomains::eGLL_Lagrange, Nektar::SpatialDomains::eGLL_Lagrange_SEM, Nektar::LibUtilities::eHexahedron, Nektar::SpatialDomains::eModified, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eModified_C, Nektar::SpatialDomains::eModifiedQuadPlus1, Nektar::SpatialDomains::eModifiedQuadPlus2, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eOrtho_C, Nektar::SpatialDomains::eOrthogonal, Nektar::LibUtilities::ePrism, Nektar::LibUtilities::ePyramid, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eSegment, Nektar::LibUtilities::eTetrahedron, and Nektar::LibUtilities::eTriangle.

Referenced by ReadExpansions().

{
LibUtilities::ShapeType shape= in->GetShapeType();
int quadoffset = 1;
switch(type)
{
case eModified:
quadoffset = 1;
break;
quadoffset = 2;
break;
quadoffset = 3;
break;
default:
break;
}
switch(type)
{
case eModified:
{
switch (shape)
{
{
const LibUtilities::PointsKey pkey(nummodes+quadoffset, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eModified_A, nummodes, pkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes+quadoffset, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eModified_A, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes+quadoffset, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eModified_A, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes+quadoffset, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eModified_A, nummodes, pkey);
returnval.push_back(bkey);
const LibUtilities::PointsKey pkey1(nummodes+quadoffset-1, LibUtilities::eGaussRadauMAlpha1Beta0);
LibUtilities::BasisKey bkey1(LibUtilities::eModified_B, nummodes, pkey1);
returnval.push_back(bkey1);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes+quadoffset, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eModified_A, nummodes, pkey);
returnval.push_back(bkey);
const LibUtilities::PointsKey pkey1(nummodes+quadoffset-1, LibUtilities::eGaussRadauMAlpha1Beta0);
LibUtilities::BasisKey bkey1(LibUtilities::eModified_B, nummodes, pkey1);
returnval.push_back(bkey1);
const LibUtilities::PointsKey pkey2(nummodes+quadoffset-1, LibUtilities::eGaussRadauMAlpha2Beta0);
LibUtilities::BasisKey bkey2(LibUtilities::eModified_C, nummodes, pkey2);
returnval.push_back(bkey2);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes+quadoffset, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eModified_A, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
const LibUtilities::PointsKey pkey1(nummodes+quadoffset-1, LibUtilities::eGaussRadauMAlpha2Beta0);
LibUtilities::BasisKey bkey1(LibUtilities::eModified_C, nummodes, pkey1);
returnval.push_back(bkey1);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes+quadoffset, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eModified_A, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
const LibUtilities::PointsKey pkey1(nummodes+quadoffset-1, LibUtilities::eGaussRadauMAlpha1Beta0);
LibUtilities::BasisKey bkey1(LibUtilities::eModified_B, nummodes, pkey1);
returnval.push_back(bkey1);
}
break;
default:
{
ASSERTL0(false,"Expansion not defined in switch for this shape");
}
break;
}
}
break;
{
switch(shape)
{
{
const LibUtilities::PointsKey pkey(nummodes+1, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eGLL_Lagrange, nummodes, pkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes+1, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eGLL_Lagrange, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
case LibUtilities::eTriangle: // define with corrects points key
// and change to Ortho on construction
{
const LibUtilities::PointsKey pkey(nummodes+1, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eGLL_Lagrange, nummodes, pkey);
returnval.push_back(bkey);
const LibUtilities::PointsKey pkey1(nummodes, LibUtilities::eGaussRadauMAlpha1Beta0);
LibUtilities::BasisKey bkey1(LibUtilities::eOrtho_B, nummodes, pkey1);
returnval.push_back(bkey1);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes+1,LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eGLL_Lagrange, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
default:
{
ASSERTL0(false, "Expansion not defined in switch for this shape");
}
break;
}
}
break;
{
switch (shape)
{
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eGaussGaussLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eGauss_Lagrange, nummodes, pkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes,LibUtilities::eGaussGaussLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eGauss_Lagrange, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes,LibUtilities::eGaussGaussLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eGauss_Lagrange, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
default:
{
ASSERTL0(false, "Expansion not defined in switch for this shape");
}
break;
}
}
break;
{
switch (shape)
{
{
const LibUtilities::PointsKey pkey(nummodes+1, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eOrtho_A, nummodes, pkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes+1, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eOrtho_A, nummodes, pkey);
returnval.push_back(bkey);
const LibUtilities::PointsKey pkey1(nummodes, LibUtilities::eGaussRadauMAlpha1Beta0);
LibUtilities::BasisKey bkey1(LibUtilities::eOrtho_B, nummodes, pkey1);
returnval.push_back(bkey1);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes+1, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eOrtho_A, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes+1, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eOrtho_A, nummodes, pkey);
returnval.push_back(bkey);
const LibUtilities::PointsKey pkey1(nummodes, LibUtilities::eGaussRadauMAlpha1Beta0);
LibUtilities::BasisKey bkey1(LibUtilities::eOrtho_B, nummodes, pkey1);
returnval.push_back(bkey1);
const LibUtilities::PointsKey pkey2(nummodes, LibUtilities::eGaussRadauMAlpha2Beta0);
LibUtilities::BasisKey bkey2(LibUtilities::eOrtho_C, nummodes, pkey2);
}
break;
default:
{
ASSERTL0(false,"Expansion not defined in switch for this shape");
}
break;
}
}
break;
{
switch (shape)
{
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eGLL_Lagrange, nummodes, pkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eGLL_Lagrange, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey(LibUtilities::eGLL_Lagrange, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
default:
{
ASSERTL0(false,"Expansion not defined in switch for this shape");
}
break;
}
}
break;
case eFourier:
{
switch (shape)
{
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eFourierEvenlySpaced);
LibUtilities::BasisKey bkey(LibUtilities::eFourier, nummodes, pkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eFourierEvenlySpaced);
LibUtilities::BasisKey bkey(LibUtilities::eFourier, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eFourierEvenlySpaced);
LibUtilities::BasisKey bkey(LibUtilities::eFourier, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
default:
{
ASSERTL0(false,"Expansion not defined in switch for this shape");
}
break;
}
}
break;
{
switch (shape)
{
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey(LibUtilities::eFourierSingleMode, nummodes, pkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey(LibUtilities::eFourierSingleMode, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey(LibUtilities::eFourierSingleMode, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
default:
{
ASSERTL0(false,"Expansion not defined in switch for this shape");
}
break;
}
}
break;
{
switch (shape)
{
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey(LibUtilities::eFourierHalfModeRe, nummodes, pkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey(LibUtilities::eFourierHalfModeRe, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey(LibUtilities::eFourierHalfModeRe, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
default:
{
ASSERTL0(false,"Expansion not defined in switch for this shape");
}
break;
}
}
break;
{
switch (shape)
{
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey(LibUtilities::eFourierHalfModeIm, nummodes, pkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey(LibUtilities::eFourierHalfModeIm, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey(LibUtilities::eFourierHalfModeIm, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
default:
{
ASSERTL0(false,"Expansion not defined in switch for this shape");
}
break;
}
}
break;
case eChebyshev:
{
switch (shape)
{
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eGaussGaussChebyshev);
LibUtilities::BasisKey bkey(LibUtilities::eChebyshev, nummodes, pkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eGaussGaussChebyshev);
LibUtilities::BasisKey bkey(LibUtilities::eChebyshev, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eGaussGaussChebyshev);
LibUtilities::BasisKey bkey(LibUtilities::eChebyshev, nummodes, pkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
returnval.push_back(bkey);
}
break;
default:
{
ASSERTL0(false,"Expansion not defined in switch for this shape");
}
break;
}
}
break;
{
switch (shape)
{
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eFourierEvenlySpaced);
LibUtilities::BasisKey bkey(LibUtilities::eFourier, nummodes, pkey);
returnval.push_back(bkey);
const LibUtilities::PointsKey pkey1(nummodes, LibUtilities::eGaussGaussChebyshev);
LibUtilities::BasisKey bkey1(LibUtilities::eChebyshev, nummodes, pkey1);
returnval.push_back(bkey1);
}
break;
default:
{
ASSERTL0(false,"Expansion not defined in switch for this shape");
}
break;
}
}
break;
{
switch (shape)
{
{
const LibUtilities::PointsKey pkey1(nummodes, LibUtilities::eGaussGaussChebyshev);
LibUtilities::BasisKey bkey1(LibUtilities::eChebyshev, nummodes, pkey1);
returnval.push_back(bkey1);
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eFourierEvenlySpaced);
LibUtilities::BasisKey bkey(LibUtilities::eFourier, nummodes, pkey);
returnval.push_back(bkey);
}
break;
default:
{
ASSERTL0(false,"Expansion not defined in switch for this shape");
}
break;
}
}
break;
{
switch (shape)
{
{
const LibUtilities::PointsKey pkey(nummodes, LibUtilities::eFourierEvenlySpaced);
LibUtilities::BasisKey bkey(LibUtilities::eFourier, nummodes, pkey);
returnval.push_back(bkey);
const LibUtilities::PointsKey pkey1(nummodes+1, LibUtilities::eGaussLobattoLegendre);
LibUtilities::BasisKey bkey1(LibUtilities::eModified_A, nummodes, pkey1);
returnval.push_back(bkey1);
}
break;
default:
{
ASSERTL0(false,"Expansion not defined in switch for this shape");
}
break;
}
}
break;
default:
{
ASSERTL0(false,"Expansion type not defined");
}
break;
}
return returnval;
}
LibUtilities::BasisKeyVector Nektar::SpatialDomains::MeshGraph::DefineBasisKeyFromExpansionTypeHomo ( GeometrySharedPtr  in,
ExpansionType  type_x,
ExpansionType  type_y,
ExpansionType  type_z,
const int  nummodes_x,
const int  nummodes_y,
const int  nummodes_z 
)

Definition at line 3167 of file MeshGraph.cpp.

References ASSERTL0, Nektar::LibUtilities::eChebyshev, Nektar::SpatialDomains::eChebyshev, Nektar::LibUtilities::eFourier, Nektar::SpatialDomains::eFourier, Nektar::LibUtilities::eFourierEvenlySpaced, Nektar::LibUtilities::eFourierHalfModeIm, Nektar::SpatialDomains::eFourierHalfModeIm, Nektar::LibUtilities::eFourierHalfModeRe, Nektar::SpatialDomains::eFourierHalfModeRe, Nektar::LibUtilities::eFourierSingleMode, Nektar::SpatialDomains::eFourierSingleMode, Nektar::LibUtilities::eFourierSingleModeSpaced, Nektar::LibUtilities::eGaussGaussChebyshev, Nektar::LibUtilities::eHexahedron, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eSegment, Nektar::LibUtilities::eTetrahedron, and Nektar::LibUtilities::eTriangle.

Referenced by ReadExpansions().

{
LibUtilities::ShapeType shape = in->GetShapeType();
switch (shape)
{
{
ASSERTL0(false,"Homogeneous expansion not defined for this shape");
}
break;
{
ASSERTL0(false,"Homogeneous expansion not defined for this shape");
}
break;
{
switch(type_x)
{
case eFourier:
{
const LibUtilities::PointsKey pkey1(nummodes_x,LibUtilities::eFourierEvenlySpaced);
LibUtilities::BasisKey bkey1(LibUtilities::eFourier,nummodes_x,pkey1);
returnval.push_back(bkey1);
}
break;
{
const LibUtilities::PointsKey pkey1(nummodes_x,LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey1(LibUtilities::eFourierSingleMode,nummodes_x,pkey1);
returnval.push_back(bkey1);
}
break;
{
const LibUtilities::PointsKey pkey1(nummodes_x,LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey1(LibUtilities::eFourierHalfModeRe,nummodes_x,pkey1);
returnval.push_back(bkey1);
}
break;
{
const LibUtilities::PointsKey pkey1(nummodes_x,LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey1(LibUtilities::eFourierHalfModeIm,nummodes_x,pkey1);
returnval.push_back(bkey1);
}
break;
case eChebyshev:
{
const LibUtilities::PointsKey pkey1(nummodes_x,LibUtilities::eGaussGaussChebyshev);
LibUtilities::BasisKey bkey1(LibUtilities::eChebyshev,nummodes_x,pkey1);
returnval.push_back(bkey1);
}
break;
default:
{
ASSERTL0(false,"Homogeneous expansion can be of Fourier or Chebyshev type only");
}
break;
}
switch(type_y)
{
case eFourier:
{
const LibUtilities::PointsKey pkey2(nummodes_y,LibUtilities::eFourierEvenlySpaced);
LibUtilities::BasisKey bkey2(LibUtilities::eFourier,nummodes_y,pkey2);
returnval.push_back(bkey2);
}
break;
{
const LibUtilities::PointsKey pkey2(nummodes_y,LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey2(LibUtilities::eFourierSingleMode,nummodes_y,pkey2);
returnval.push_back(bkey2);
}
break;
{
const LibUtilities::PointsKey pkey2(nummodes_y,LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey2(LibUtilities::eFourierHalfModeRe,nummodes_y,pkey2);
returnval.push_back(bkey2);
}
break;
{
const LibUtilities::PointsKey pkey2(nummodes_y,LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey2(LibUtilities::eFourierHalfModeIm,nummodes_y,pkey2);
returnval.push_back(bkey2);
}
break;
case eChebyshev:
{
const LibUtilities::PointsKey pkey2(nummodes_y,LibUtilities::eGaussGaussChebyshev);
LibUtilities::BasisKey bkey2(LibUtilities::eChebyshev,nummodes_y,pkey2);
returnval.push_back(bkey2);
}
break;
default:
{
ASSERTL0(false,"Homogeneous expansion can be of Fourier or Chebyshev type only");
}
break;
}
switch(type_z)
{
case eFourier:
{
const LibUtilities::PointsKey pkey3(nummodes_z,LibUtilities::eFourierEvenlySpaced);
LibUtilities::BasisKey bkey3(LibUtilities::eFourier,nummodes_z,pkey3);
returnval.push_back(bkey3);
}
break;
{
const LibUtilities::PointsKey pkey3(nummodes_z,LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey3(LibUtilities::eFourierSingleMode,nummodes_z,pkey3);
returnval.push_back(bkey3);
}
break;
{
const LibUtilities::PointsKey pkey3(nummodes_z,LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey3(LibUtilities::eFourierHalfModeRe,nummodes_z,pkey3);
returnval.push_back(bkey3);
}
break;
{
const LibUtilities::PointsKey pkey3(nummodes_z,LibUtilities::eFourierSingleModeSpaced);
LibUtilities::BasisKey bkey3(LibUtilities::eFourierHalfModeIm,nummodes_z,pkey3);
returnval.push_back(bkey3);
}
break;
case eChebyshev:
{
const LibUtilities::PointsKey pkey3(nummodes_z,LibUtilities::eGaussGaussChebyshev);
LibUtilities::BasisKey bkey3(LibUtilities::eChebyshev,nummodes_z,pkey3);
returnval.push_back(bkey3);
}
break;
default:
{
ASSERTL0(false,"Homogeneous expansion can be of Fourier or Chebyshev type only");
}
break;
}
}
break;
{
ASSERTL0(false,"Homogeneous expansion not defined for this shape");
}
break;
{
ASSERTL0(false,"Homogeneous expansion not defined for this shape");
}
break;
default:
ASSERTL0(false,"Expansion not defined in switch for this shape");
break;
}
return returnval;
}
const std::map< int, boost::shared_ptr< PyrGeom > > & Nektar::SpatialDomains::MeshGraph::GetAllElementsOfType ( ) const
inline

Convenience method for ElVis.

Definition at line 567 of file MeshGraph.h.

References GetAllHexGeoms().

{
return GetAllHexGeoms();
}
const HexGeomMap& Nektar::SpatialDomains::MeshGraph::GetAllHexGeoms ( ) const
inline

Definition at line 378 of file MeshGraph.h.

References m_hexGeoms.

Referenced by GetAllElementsOfType().

{ return m_hexGeoms; }
const PrismGeomMap& Nektar::SpatialDomains::MeshGraph::GetAllPrismGeoms ( ) const
inline

Definition at line 377 of file MeshGraph.h.

References m_prismGeoms.

{ return m_prismGeoms; }
const PyrGeomMap& Nektar::SpatialDomains::MeshGraph::GetAllPyrGeoms ( ) const
inline

Definition at line 376 of file MeshGraph.h.

References m_pyrGeoms.

{ return m_pyrGeoms; }
const QuadGeomMap& Nektar::SpatialDomains::MeshGraph::GetAllQuadGeoms ( ) const
inline

Definition at line 374 of file MeshGraph.h.

References m_quadGeoms.

{ return m_quadGeoms; }
const SegGeomMap& Nektar::SpatialDomains::MeshGraph::GetAllSegGeoms ( ) const
inline

Definition at line 372 of file MeshGraph.h.

References m_segGeoms.

{ return m_segGeoms; }
const TetGeomMap& Nektar::SpatialDomains::MeshGraph::GetAllTetGeoms ( ) const
inline

Definition at line 375 of file MeshGraph.h.

References m_tetGeoms.

{ return m_tetGeoms; }
const TriGeomMap& Nektar::SpatialDomains::MeshGraph::GetAllTriGeoms ( ) const
inline

Definition at line 373 of file MeshGraph.h.

References m_triGeoms.

{ return m_triGeoms; }
Composite Nektar::SpatialDomains::MeshGraph::GetComposite ( int  whichComposite) const
inline

Definition at line 441 of file MeshGraph.h.

References ASSERTL0, and m_meshComposites.

Referenced by GetCompositeList(), and Nektar::SpatialDomains::Domain::Read().

{
Composite returnval;
ASSERTL0(m_meshComposites.find(whichComposite) != m_meshComposites.end(),
"Composite not found.");
return m_meshComposites.find(whichComposite)->second;
}
GeometrySharedPtr Nektar::SpatialDomains::MeshGraph::GetCompositeItem ( int  whichComposite,
int  whichItem 
)

Definition at line 1616 of file MeshGraph.cpp.

References ErrorUtil::efatal, m_meshComposites, and NEKERROR.

{
GeometrySharedPtr returnval;
bool error = false;
if (whichComposite >= 0 && whichComposite < int(m_meshComposites.size()))
{
if (whichItem >= 0 && whichItem < int(m_meshComposites[whichComposite]->size()))
{
returnval = m_meshComposites[whichComposite]->at(whichItem);
}
else
{
error = true;
}
}
else
{
error = true;
}
if (error)
{
std::ostringstream errStream;
errStream << "Unable to access composite item [" << whichComposite << "][" << whichItem << "].";
std::string testStr = errStream.str();
NEKERROR(ErrorUtil::efatal, testStr.c_str());
}
return returnval;
}
void Nektar::SpatialDomains::MeshGraph::GetCompositeList ( const std::string &  compositeStr,
CompositeMap compositeVector 
) const

Definition at line 1654 of file MeshGraph.cpp.

References ASSERTL0, ErrorUtil::ewarning, Nektar::StdRegions::find(), Nektar::ParseUtils::GenerateSeqVector(), GetComposite(), Nektar::iterator, m_meshComposites, m_meshPartitioned, and NEKERROR.

Referenced by ReadDomain(), and ReadExpansions().

{
// Parse the composites into a list.
typedef vector<unsigned int> SeqVector;
SeqVector seqVector;
bool parseGood = ParseUtils::GenerateSeqVector(compositeStr.c_str(), seqVector);
ASSERTL0(parseGood && !seqVector.empty(), (std::string("Unable to read composite index range: ") + compositeStr).c_str());
SeqVector addedVector; // Vector of those composites already added to compositeVector;
for (SeqVector::iterator iter = seqVector.begin(); iter != seqVector.end(); ++iter)
{
// Only add a new one if it does not already exist in vector.
// Can't go back and delete with a vector, so prevent it from
// being added in the first place.
if (std::find(addedVector.begin(), addedVector.end(), *iter) == addedVector.end())
{
// If the composite listed is not found and we are working
// on a partitioned mesh, silently ignore it.
if (m_meshComposites.find(*iter) == m_meshComposites.end()
{
continue;
}
addedVector.push_back(*iter);
Composite composite = GetComposite(*iter);
if (composite)
{
compositeVector[*iter] = composite;
}
else
{
char str[64];
::sprintf(str, "%d", *iter);
NEKERROR(ErrorUtil::ewarning, (std::string("Undefined composite: ") + str).c_str());
}
}
}
}
const CompositeMap & Nektar::SpatialDomains::MeshGraph::GetComposites ( ) const
inline

Definition at line 453 of file MeshGraph.h.

References m_meshComposites.

{
}
const CurveVector& Nektar::SpatialDomains::MeshGraph::GetCurvedEdges ( ) const
inline

Definition at line 368 of file MeshGraph.h.

References m_curvedEdges.

{ return m_curvedEdges; }
const CurveVector& Nektar::SpatialDomains::MeshGraph::GetCurvedFaces ( ) const
inline

Definition at line 370 of file MeshGraph.h.

References m_curvedFaces.

{ return m_curvedFaces; }
const std::vector< CompositeMap > & Nektar::SpatialDomains::MeshGraph::GetDomain ( void  ) const
inline

Definition at line 462 of file MeshGraph.h.

References m_domain.

{
return m_domain;
}
const CompositeMap & Nektar::SpatialDomains::MeshGraph::GetDomain ( int  domain) const
inline

Definition at line 470 of file MeshGraph.h.

References ASSERTL1, and m_domain.

{
ASSERTL1(domain < m_domain.size(),"Request for domain which does not exist");
return m_domain[domain];
}
SegGeomSharedPtr Nektar::SpatialDomains::MeshGraph::GetEdge ( unsigned int  id)
inline

Definition at line 357 of file MeshGraph.h.

References m_segGeoms.

{ return m_segGeoms[id]; }
ExpansionShPtr Nektar::SpatialDomains::MeshGraph::GetExpansion ( GeometrySharedPtr  geom,
const std::string  variable = "DefaultVar" 
)

Definition at line 1729 of file MeshGraph.cpp.

References m_expansionMapShPtrMap.

Referenced by Nektar::SpatialDomains::MeshGraph2D::GetEdgeBasisKey(), and Nektar::SpatialDomains::MeshGraph3D::GetFaceBasisKey().

{
ExpansionShPtr returnval;
ExpansionMapShPtr expansionMap = m_expansionMapShPtrMap.find(variable)->second;
for (iter = expansionMap->begin(); iter!=expansionMap->end(); ++iter)
{
if ((iter->second)->m_geomShPtr == geom)
{
returnval = iter->second;
break;
}
}
return returnval;
}
const ExpansionMap & Nektar::SpatialDomains::MeshGraph::GetExpansions ( )
inline

Definition at line 480 of file MeshGraph.h.

{
std::string defstr = "DefaultVar";
return GetExpansions(defstr);
}
const ExpansionMap & Nektar::SpatialDomains::MeshGraph::GetExpansions ( const std::string  variable)

Definition at line 1702 of file MeshGraph.cpp.

References ErrorUtil::efatal, ErrorUtil::ewarning, m_expansionMapShPtrMap, and NEKERROR.

{
ExpansionMapShPtr returnval;
if(m_expansionMapShPtrMap.count(variable))
{
returnval = m_expansionMapShPtrMap.find(variable)->second;
}
else
{
if(m_expansionMapShPtrMap.count("DefaultVar") == 0)
{
NEKERROR(ErrorUtil::efatal, (std::string("Unable to find expansion vector definition for field: ")+variable).c_str());
}
returnval = m_expansionMapShPtrMap.find("DefaultVar")->second;
m_expansionMapShPtrMap[variable] = returnval;
NEKERROR(ErrorUtil::ewarning, (std::string("Using Default variable expansion definition for field: ")+variable).c_str());
}
return *returnval;
}
const std::string Nektar::SpatialDomains::MeshGraph::GetGeomInfo ( std::string  parameter)
inline

Definition at line 532 of file MeshGraph.h.

References ASSERTL1, and m_geomInfo.

{
ASSERTL1(m_geomInfo.find(parameter) != m_geomInfo.end(),
"Parameter " + parameter + " does not exist.");
return m_geomInfo[parameter];
}
int Nektar::SpatialDomains::MeshGraph::GetMeshDimension ( void  ) const
inline

Dimension of the mesh (can be a 1D curve in 3D space).

Definition at line 423 of file MeshGraph.h.

References m_meshDimension.

Referenced by Read().

{
}
int Nektar::SpatialDomains::MeshGraph::GetNvertices ( ) const
inline

Definition at line 543 of file MeshGraph.h.

References m_vertSet.

{
return int(m_vertSet.size());
}
int Nektar::SpatialDomains::MeshGraph::GetSpaceDimension ( void  ) const
inline

Dimension of the space (can be a 1D curve in 3D space).

Definition at line 432 of file MeshGraph.h.

References m_spaceDimension.

Referenced by Nektar::SpatialDomains::MeshGraph1D::GetCoordim(), Nektar::SpatialDomains::MeshGraph3D::GetCoordim(), and Nektar::SpatialDomains::MeshGraph2D::GetCoordim().

{
}
PointGeomSharedPtr Nektar::SpatialDomains::MeshGraph::GetVertex ( int  id)
inline

Definition at line 552 of file MeshGraph.h.

References ASSERTL0, Nektar::iterator, and m_vertSet.

Referenced by Nektar::SpatialDomains::MeshGraph3D::ReadEdges(), Nektar::SpatialDomains::MeshGraph2D::ReadEdges(), and Nektar::SpatialDomains::MeshGraph1D::ReadElements().

{
PointGeomSharedPtr returnval;
ASSERTL0(x != m_vertSet.end(),
"Vertex " + boost::lexical_cast<string>(id)
+ " not found.");
return x->second;
}
boost::shared_ptr< MeshGraph > Nektar::SpatialDomains::MeshGraph::Read ( const LibUtilities::SessionReaderSharedPtr pSession,
DomainRangeShPtr rng = NullDomainRangeShPtr 
)
static

Definition at line 118 of file MeshGraph.cpp.

References ASSERTL1, ErrorUtil::efatal, and NEKERROR.

Referenced by Diffusion::Diffusion(), Nektar::VortexWaveInteraction::FileRelaxation(), main(), Nektar::Utilities::InputNekpp::Process(), Nektar::Utilities::InputXml::Process(), Nektar::Utilities::ProcessInterpPoints::Process(), Nektar::Utilities::ProcessInterpField::Process(), Nektar::SolverUtils::EquationSystem::v_InitObject(), and Nektar::SolverUtils::FilterModalEnergy::v_Update().

{
boost::shared_ptr<MeshGraph> returnval;
// read the geometry tag to get the dimension
TiXmlElement* geometry_tag = pSession->GetElement("NEKTAR/GEOMETRY");
TiXmlAttribute *attr = geometry_tag->FirstAttribute();
int meshDim = 0;
while (attr)
{
std::string attrName(attr->Name());
if (attrName == "DIM")
{
int err = attr->QueryIntValue(&meshDim);
ASSERTL1(err==TIXML_SUCCESS, "Unable to read mesh dimension.");
break;
}
else
{
std::string errstr("Unknown attribute: ");
errstr += attrName;
ASSERTL1(false, errstr.c_str());
}
// Get the next attribute.
attr = attr->Next();
}
// instantiate the dimension-specific meshgraph classes
switch(meshDim)
{
case 1:
returnval = MemoryManager<MeshGraph1D>::AllocateSharedPtr(pSession,rng);
break;
case 2:
returnval = MemoryManager<MeshGraph2D>::AllocateSharedPtr(pSession,rng);
break;
case 3:
returnval = MemoryManager<MeshGraph3D>::AllocateSharedPtr(pSession,rng);
break;
default:
std::string err = "Invalid mesh dimension: ";
std::stringstream strstrm;
strstrm << meshDim;
err += strstrm.str();
NEKERROR(ErrorUtil::efatal, err.c_str());
}
return returnval;
}
boost::shared_ptr< MeshGraph > Nektar::SpatialDomains::MeshGraph::Read ( const std::string &  infilename,
bool  pReadExpansions = true 
)
static
Todo:
Remove updated routine

Definition at line 179 of file MeshGraph.cpp.

References ErrorUtil::efatal, GetMeshDimension(), NEKERROR, and ReadGeometry().

{
boost::shared_ptr<MeshGraph> returnval;
MeshGraph mesh;
mesh.ReadGeometry(infilename);
int meshDim = mesh.GetMeshDimension();
switch(meshDim)
{
case 1:
returnval = MemoryManager<MeshGraph1D>::AllocateSharedPtr();
break;
case 2:
returnval = MemoryManager<MeshGraph2D>::AllocateSharedPtr();
break;
case 3:
returnval = MemoryManager<MeshGraph3D>::AllocateSharedPtr();
break;
default:
std::string err = "Invalid mesh dimension: ";
std::stringstream strstrm;
strstrm << meshDim;
err += strstrm.str();
NEKERROR(ErrorUtil::efatal, err.c_str());
}
if (returnval)
{
returnval->ReadGeometry(infilename);
returnval->ReadGeometryInfo(infilename);
if (pReadExpansions)
{
returnval->ReadExpansions(infilename);
}
}
return returnval;
}
void Nektar::SpatialDomains::MeshGraph::ReadCurves ( TiXmlDocument &  doc)

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

Look for elements in CURVE block.

All curves are of the form: "<? ID="#" TYPE="GLL OR other points type" NUMPOINTS="#"> ... </?>", with ? being an element type (either E or F).

These should be ordered.

Read id attribute.

Read edge id attribute.

Read text edgelement description.

Parse out the element components corresponding to type of element.

Read id attribute.

Read face id attribute.

Read text face element description.

Parse out the element components corresponding to type of element.

Definition at line 1102 of file MeshGraph.cpp.

References ASSERTL0, Nektar::LibUtilities::AnalyticExpressionEvaluator::DefineFunction(), ErrorUtil::efatal, Nektar::LibUtilities::AnalyticExpressionEvaluator::Evaluate(), Nektar::StdRegions::find(), Nektar::LibUtilities::kPointsTypeStr, m_curvedEdges, m_curvedFaces, m_meshDimension, NEKERROR, and Nektar::LibUtilities::SIZE_PointsType.

Referenced by ReadCurves(), Nektar::SpatialDomains::MeshGraph3D::ReadGeometry(), and Nektar::SpatialDomains::MeshGraph2D::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;
// check to see if any scaling parameters are in
// attributes and determine these values
TiXmlElement* element = mesh->FirstChildElement("VERTEX");
ASSERTL0(element, "Unable to find mesh VERTEX tag in file.");
NekDouble xscale,yscale,zscale;
LibUtilities::AnalyticExpressionEvaluator expEvaluator;
const char *xscal = element->Attribute("XSCALE");
if(!xscal)
{
xscale = 1.0;
}
else
{
std::string xscalstr = xscal;
int expr_id = expEvaluator.DefineFunction("",xscalstr);
xscale = expEvaluator.Evaluate(expr_id);
}
const char *yscal = element->Attribute("YSCALE");
if(!yscal)
{
yscale = 1.0;
}
else
{
std::string yscalstr = yscal;
int expr_id = expEvaluator.DefineFunction("",yscalstr);
yscale = expEvaluator.Evaluate(expr_id);
}
const char *zscal = element->Attribute("ZSCALE");
if(!zscal)
{
zscale = 1.0;
}
else
{
std::string zscalstr = zscal;
int expr_id = expEvaluator.DefineFunction("",zscalstr);
zscale = expEvaluator.Evaluate(expr_id);
}
int err;
/// Look for elements in CURVE block.
field = mesh->FirstChildElement("CURVED");
if(!field) //return if no curved entities
{
return;
}
/// All curves are of the form: "<? ID="#" TYPE="GLL OR other
/// points type" NUMPOINTS="#"> ... </?>", with ? being an
/// element type (either E or F).
TiXmlElement *edgelement = field->FirstChildElement("E");
int edgeindx, edgeid;
int nextEdgeNumber = -1;
while(edgelement)
{
/// These should be ordered.
nextEdgeNumber++;
std::string edge(edgelement->ValueStr());
ASSERTL0(edge == "E", (std::string("Unknown 3D curve type:") + edge).c_str());
/// Read id attribute.
err = edgelement->QueryIntAttribute("ID", &edgeindx);
ASSERTL0(err == TIXML_SUCCESS, "Unable to read curve attribute ID.");
/// Read edge id attribute.
err = edgelement->QueryIntAttribute("EDGEID", &edgeid);
ASSERTL0(err == TIXML_SUCCESS, "Unable to read curve attribute EDGEID.");
/// Read text edgelement description.
std::string elementStr;
TiXmlNode* elementChild = edgelement->FirstChild();
while(elementChild)
{
// Accumulate all non-comment element data
if (elementChild->Type() == TiXmlNode::TEXT)
{
elementStr += elementChild->ToText()->ValueStr();
elementStr += " ";
}
elementChild = elementChild->NextSibling();
}
ASSERTL0(!elementStr.empty(), "Unable to read curve description body.");
/// Parse out the element components corresponding to type of element.
if (edge == "E")
{
int numPts=0;
// Determine the points type
std::string typeStr = edgelement->Attribute("TYPE");
ASSERTL0(!typeStr.empty(), "TYPE must be specified in " "points definition");
const std::string* begStr = LibUtilities::kPointsTypeStr;
const std::string* ptsStr = std::find(begStr, endStr, typeStr);
ASSERTL0(ptsStr != endStr, "Invalid points type.");
type = (LibUtilities::PointsType)(ptsStr - begStr);
//Determine the number of points
err = edgelement->QueryIntAttribute("NUMPOINTS", &numPts);
ASSERTL0(err == TIXML_SUCCESS, "Unable to read curve attribute NUMPOINTS.");
CurveSharedPtr curve(MemoryManager<Curve>::AllocateSharedPtr(edgeid, type));
// Read points (x, y, z)
NekDouble xval, yval, zval;
std::istringstream elementDataStrm(elementStr.c_str());
try
{
while(!elementDataStrm.fail())
{
elementDataStrm >> xval >> yval >> zval;
xval *= xscale;
yval *= yscale;
zval *= zscale;
// Need to check it here because we may not be
// good after the read indicating that there
// was nothing to read.
if (!elementDataStrm.fail())
{
PointGeomSharedPtr vert(MemoryManager<PointGeom>::AllocateSharedPtr(m_meshDimension, edgeindx, xval, yval, zval));
curve->m_points.push_back(vert);
}
}
}
catch(...)
{
(std::string("Unable to read curve data for EDGE: ") + elementStr).c_str());
}
ASSERTL0(curve->m_points.size() == numPts,"Number of points specificed by attribute NUMPOINTS is different from number of points in list");
m_curvedEdges.push_back(curve);
edgelement = edgelement->NextSiblingElement("E");
} // end if-loop
} // end while-loop
TiXmlElement *facelement = field->FirstChildElement("F");
int faceindx, faceid;
while(facelement)
{
std::string face(facelement->ValueStr());
ASSERTL0(face == "F", (std::string("Unknown 3D curve type: ") + face).c_str());
/// Read id attribute.
err = facelement->QueryIntAttribute("ID", &faceindx);
ASSERTL0(err == TIXML_SUCCESS, "Unable to read curve attribute ID.");
/// Read face id attribute.
err = facelement->QueryIntAttribute("FACEID", &faceid);
ASSERTL0(err == TIXML_SUCCESS, "Unable to read curve attribute FACEID.");
/// Read text face element description.
std::string elementStr;
TiXmlNode* elementChild = facelement->FirstChild();
while(elementChild)
{
// Accumulate all non-comment element data
if (elementChild->Type() == TiXmlNode::TEXT)
{
elementStr += elementChild->ToText()->ValueStr();
elementStr += " ";
}
elementChild = elementChild->NextSibling();
}
ASSERTL0(!elementStr.empty(), "Unable to read curve description body.");
/// Parse out the element components corresponding to type of element.
if(face == "F")
{
std::string typeStr = facelement->Attribute("TYPE");
ASSERTL0(!typeStr.empty(), "TYPE must be specified in " "points definition");
const std::string* begStr = LibUtilities::kPointsTypeStr;
const std::string* ptsStr = std::find(begStr, endStr, typeStr);
ASSERTL0(ptsStr != endStr, "Invalid points type.");
type = (LibUtilities::PointsType)(ptsStr - begStr);
std::string numptsStr = facelement->Attribute("NUMPOINTS");
ASSERTL0(!numptsStr.empty(), "NUMPOINTS must be specified in points definition");
int numPts=0;
std::stringstream s;
s << numptsStr;
s >> numPts;
CurveSharedPtr curve(MemoryManager<Curve>::AllocateSharedPtr(faceid, type));
ASSERTL0(numPts >= 3, "NUMPOINTS for face must be greater than 2");
if(numPts == 3)
{
ASSERTL0(ptsStr != endStr, "Invalid points type.");
}
// Read points (x, y, z)
NekDouble xval, yval, zval;
std::istringstream elementDataStrm(elementStr.c_str());
try
{
while(!elementDataStrm.fail())
{
elementDataStrm >> xval >> yval >> zval;
// Need to check it here because we may not be good after the read
// indicating that there was nothing to read.
if (!elementDataStrm.fail())
{
PointGeomSharedPtr vert(MemoryManager<PointGeom>::AllocateSharedPtr(m_meshDimension, faceindx, xval, yval, zval));
curve->m_points.push_back(vert);
}
}
}
catch(...)
{
(std::string("Unable to read curve data for FACE: ")
+ elementStr).c_str());
}
m_curvedFaces.push_back(curve);
facelement = facelement->NextSiblingElement("F");
} // end if-loop
} // end while-loop
} // end of ReadCurves()
void Nektar::SpatialDomains::MeshGraph::ReadCurves ( std::string &  infilename)

Definition at line 1366 of file MeshGraph.cpp.

References ASSERTL0, and ReadCurves().

{
TiXmlDocument doc(infilename);
bool loadOkay = doc.LoadFile();
std::stringstream errstr;
errstr << "Unable to load file: " << infilename << std::endl;
errstr << "Reason: " << doc.ErrorDesc() << std::endl;
errstr << "Position: Line " << doc.ErrorRow() << ", Column " << doc.ErrorCol() << std::endl;
ASSERTL0(loadOkay, errstr.str());
ReadCurves(doc);
}
void Nektar::SpatialDomains::MeshGraph::ReadDomain ( TiXmlDocument &  doc)

Look for data in DOMAIN block.

Elements are of the form: "<D ID = "N"> ... </D>". Read the ID field first.

Keep looking

Definition at line 1008 of file MeshGraph.cpp.

References ASSERTL0, GetCompositeList(), and m_domain.

Referenced by Nektar::SpatialDomains::MeshGraph1D::ReadGeometry(), Nektar::SpatialDomains::MeshGraph3D::ReadGeometry(), and Nektar::SpatialDomains::MeshGraph2D::ReadGeometry().

{
TiXmlHandle docHandle(&doc);
TiXmlElement* mesh = docHandle.FirstChildElement("NEKTAR").FirstChildElement("GEOMETRY").Element();
TiXmlElement* domain = NULL;
ASSERTL0(mesh, "Unable to find GEOMETRY tag in file.");
/// Look for data in DOMAIN block.
domain = mesh->FirstChildElement("DOMAIN");
ASSERTL0(domain, "Unable to find DOMAIN tag in file.");
/// Elements are of the form: "<D ID = "N"> ... </D>".
/// Read the ID field first.
TiXmlElement *multidomains = domain->FirstChildElement("D");
if(multidomains)
{
int nextDomainNumber = 0;
while (multidomains)
{
int indx;
int err = multidomains->QueryIntAttribute("ID", &indx);
ASSERTL0(err == TIXML_SUCCESS,
"Unable to read attribute ID in Domain.");
TiXmlNode* elementChild = multidomains->FirstChild();
while(elementChild && elementChild->Type() != TiXmlNode::TEXT)
{
elementChild = elementChild->NextSibling();
}
ASSERTL0(elementChild, "Unable to read DOMAIN body.");
std::string elementStr = elementChild->ToText()->ValueStr();
elementStr = elementStr.substr(elementStr.find_first_not_of(" "));
std::string::size_type indxBeg = elementStr.find_first_of('[') + 1;
std::string::size_type indxEnd = elementStr.find_last_of(']') - 1;
std::string indxStr = elementStr.substr(indxBeg, indxEnd - indxBeg + 1);
ASSERTL0(!indxStr.empty(), "Unable to read domain's composite index (index missing?).");
// Read the domain composites.
// Parse the composites into a list.
CompositeMap unrollDomain;
GetCompositeList(indxStr, unrollDomain);
m_domain.push_back(unrollDomain);
ASSERTL0(!m_domain[nextDomainNumber++].empty(), (std::string("Unable to obtain domain's referenced composite: ") + indxStr).c_str());
/// Keep looking
multidomains = multidomains->NextSiblingElement("D");
}
}
else // previous definition of just one composite
{
// find the non comment portion of the body.
TiXmlNode* elementChild = domain->FirstChild();
while(elementChild && elementChild->Type() != TiXmlNode::TEXT)
{
elementChild = elementChild->NextSibling();
}
ASSERTL0(elementChild, "Unable to read DOMAIN body.");
std::string elementStr = elementChild->ToText()->ValueStr();
elementStr = elementStr.substr(elementStr.find_first_not_of(" "));
std::string::size_type indxBeg = elementStr.find_first_of('[') + 1;
std::string::size_type indxEnd = elementStr.find_last_of(']') - 1;
std::string indxStr = elementStr.substr(indxBeg, indxEnd - indxBeg + 1);
ASSERTL0(!indxStr.empty(), "Unable to read domain's composite index (index missing?).");
// Read the domain composites.
// Parse the composites into a list.
CompositeMap fullDomain;
GetCompositeList(indxStr, fullDomain);
m_domain.push_back(fullDomain);
ASSERTL0(!m_domain[0].empty(), (std::string("Unable to obtain domain's referenced composite: ") + indxStr).c_str());
}
}
void Nektar::SpatialDomains::MeshGraph::ReadExpansions ( const std::string &  infilename)

Read the expansions given the XML file path.

Definition at line 524 of file MeshGraph.cpp.

References ASSERTL0.

Referenced by Nektar::SpatialDomains::MeshGraph1D::MeshGraph1D(), Nektar::SpatialDomains::MeshGraph2D::MeshGraph2D(), and Nektar::SpatialDomains::MeshGraph3D::MeshGraph3D().

{
TiXmlDocument doc(infilename);
bool loadOkay = doc.LoadFile();
std::stringstream errstr;
errstr << "Unable to load file: " << infilename << std::endl;
errstr << "Reason: " << doc.ErrorDesc() << std::endl;
errstr << "Position: Line " << doc.ErrorRow() << ", Column " << doc.ErrorCol() << std::endl;
ASSERTL0(loadOkay, errstr.str());
}
void Nektar::SpatialDomains::MeshGraph::ReadExpansions ( TiXmlDocument &  doc)

Read the expansions given the XML document reference.

Expansiontypes will contain composite, nummodes, and expansiontype (eModified, or eOrthogonal) Or a full list of data of basistype, nummodes, pointstype, numpoints;

Expansiontypes may also contain a list of fields that this expansion relates to. If this does not exist the variable is only set to "DefaultVar".

Mandatory components...optional are to follow later.

Todo:
solvers break the pattern 'instantiate Session -> instantiate MeshGraph' and parse command line arguments by themselves; one needs to unify command line arguments handling. Solvers tend to call MeshGraph::Read statically -> m_session is not defined -> no info about command line arguments presented ASSERTL0(m_session != 0, "One needs to instantiate SessionReader first");

Mandatory components...optional are to follow later.

Definition at line 542 of file MeshGraph.cpp.

References ASSERTL0, Nektar::LibUtilities::BasisTypeMap, DefineBasisKeyFromExpansionType(), DefineBasisKeyFromExpansionTypeHomo(), Nektar::SpatialDomains::eExpansionTypeSize, Nektar::SpatialDomains::eNoExpansionType, Nektar::LibUtilities::Equation::Evaluate(), Nektar::StdRegions::find(), Nektar::ParseUtils::GenerateOrderedStringVector(), Nektar::ParseUtils::GenerateOrderedVector(), GetCompositeList(), Nektar::LibUtilities::FieldIO::ImportFieldDefs(), Nektar::SpatialDomains::kExpansionTypeStr, Nektar::LibUtilities::kPointsTypeStr, m_expansionMapShPtrMap, m_session, SetExpansions(), SetUpExpansionMap(), Nektar::LibUtilities::SIZE_BasisType, and Nektar::LibUtilities::SIZE_PointsType.

{
TiXmlElement *master = doc.FirstChildElement("NEKTAR");
ASSERTL0(master, "Unable to find NEKTAR tag in file.");
// Find the Expansions tag
TiXmlElement *expansionTypes = master->FirstChildElement("EXPANSIONS");
ASSERTL0(expansionTypes, "Unable to find EXPANSIONS tag in file.");
if(expansionTypes)
{
// Find the Expansion type
TiXmlElement *expansion = expansionTypes->FirstChildElement();
std::string expType = expansion->Value();
if(expType == "E")
{
int i;
ExpansionMapShPtr expansionMap;
/// Expansiontypes will contain composite,
/// nummodes, and expansiontype (eModified, or
/// eOrthogonal) Or a full list of data of
/// basistype, nummodes, pointstype, numpoints;
/// Expansiontypes may also contain a list of
/// fields that this expansion relates to. If this
/// does not exist the variable is only set to
/// "DefaultVar".
while (expansion)
{
const char *fStr = expansion->Attribute("FIELDS");
std::vector<std::string> fieldStrings;
if(fStr) // extract other fields.
{
std::string fieldStr = fStr;
bool valid = ParseUtils::GenerateOrderedStringVector(fieldStr.c_str(),fieldStrings);
ASSERTL0(valid,"Unable to correctly parse the field string in ExpansionTypes.");
}
// check to see if m_expasionVectorShPtrMap has
// already been intiailised and if not intiailse
// vector.
if(m_expansionMapShPtrMap.count("DefaultVar") == 0) // no previous definitions
{
expansionMap = SetUpExpansionMap();
m_expansionMapShPtrMap["DefaultVar"] = expansionMap;
// make sure all fields in this search point
// to same expansion vector;
for(i = 0; i < fieldStrings.size(); ++i)
{
m_expansionMapShPtrMap[fieldStrings[i]] = expansionMap;
}
}
else // default variable is defined
{
if(fieldStrings.size()) // fields are defined
{
//see if field exists
if(m_expansionMapShPtrMap.count(fieldStrings[0]))
{
expansionMap = m_expansionMapShPtrMap.find(fieldStrings[0])->second;
}
else
{
expansionMap = SetUpExpansionMap();
// make sure all fields in this search point
// to same expansion vector;
for(i = 0; i < fieldStrings.size(); ++i)
{
if(m_expansionMapShPtrMap.count(fieldStrings[i]) == 0)
{
m_expansionMapShPtrMap[fieldStrings[i]] = expansionMap;
}
else
{
ASSERTL0(false,"Expansion vector for this field is already setup");
}
}
}
}
else // use default variable list
{
expansionMap = m_expansionMapShPtrMap.find("DefaultVar")->second;
}
}
/// Mandatory components...optional are to follow later.
std::string compositeStr = expansion->Attribute("COMPOSITE");
ASSERTL0(compositeStr.length() > 3, "COMPOSITE must be specified in expansion definition");
int beg = compositeStr.find_first_of("[");
int end = compositeStr.find_first_of("]");
std::string compositeListStr = compositeStr.substr(beg+1,end-beg-1);
CompositeMap compositeVector;
GetCompositeList(compositeListStr, compositeVector);
bool useExpansionType = false;
ExpansionType expansion_type;
int num_modes;
const char * tStr = expansion->Attribute("TYPE");
if(tStr) // use type string to define expansion
{
std::string typeStr = tStr;
const std::string* begStr = kExpansionTypeStr;
const std::string* endStr = kExpansionTypeStr+eExpansionTypeSize;
const std::string* expStr = std::find(begStr, endStr, typeStr);
ASSERTL0(expStr != endStr, "Invalid expansion type.");
expansion_type = (ExpansionType)(expStr - begStr);
/// \todo solvers break the pattern 'instantiate Session -> instantiate MeshGraph'
/// and parse command line arguments by themselves; one needs to unify command
/// line arguments handling.
/// Solvers tend to call MeshGraph::Read statically -> m_session
/// is not defined -> no info about command line arguments presented
/// ASSERTL0(m_session != 0, "One needs to instantiate SessionReader first");
const char *nStr = expansion->Attribute("NUMMODES");
ASSERTL0(nStr,"NUMMODES was not defined in EXPANSION section of input");
std::string nummodesStr = nStr;
// ASSERTL0(m_session,"Session should be defined to evaluate nummodes ");
if (m_session)
{
LibUtilities::Equation nummodesEqn(m_session, nummodesStr);
num_modes = (int) nummodesEqn.Evaluate();
}
else
{
num_modes = boost::lexical_cast<int>(nummodesStr);
}
useExpansionType = true;
}
else // assume expansion is defined individually
{
// Extract the attributes.
const char *bTypeStr = expansion->Attribute("BASISTYPE");
ASSERTL0(bTypeStr,"TYPE or BASISTYPE was not defined in EXPANSION section of input");
std::string basisTypeStr = bTypeStr;
// interpret the basis type string.
std::vector<std::string> basisStrings;
std::vector<LibUtilities::BasisType> basis;
bool valid = ParseUtils::GenerateOrderedStringVector(basisTypeStr.c_str(), basisStrings);
ASSERTL0(valid, "Unable to correctly parse the basis types.");
for (vector<std::string>::size_type i = 0; i < basisStrings.size(); i++)
{
valid = false;
for (unsigned int j = 0; j < LibUtilities::SIZE_BasisType; j++)
{
if (LibUtilities::BasisTypeMap[j] == basisStrings[i])
{
basis.push_back((LibUtilities::BasisType) j);
valid = true;
break;
}
}
ASSERTL0(valid, std::string("Unable to correctly parse the basis type: ").append(basisStrings[i]).c_str());
}
const char *nModesStr = expansion->Attribute("NUMMODES");
ASSERTL0(nModesStr,"NUMMODES was not defined in EXPANSION section of input");
std::string numModesStr = nModesStr;
std::vector<unsigned int> numModes;
valid = ParseUtils::GenerateOrderedVector(numModesStr.c_str(), numModes);
ASSERTL0(valid, "Unable to correctly parse the number of modes.");
ASSERTL0(numModes.size() == basis.size(),"information for num modes does not match the number of basis");
const char *pTypeStr = expansion->Attribute("POINTSTYPE");
ASSERTL0(pTypeStr,"POINTSTYPE was not defined in EXPANSION section of input");
std::string pointsTypeStr = pTypeStr;
// interpret the points type string.
std::vector<std::string> pointsStrings;
std::vector<LibUtilities::PointsType> points;
valid = ParseUtils::GenerateOrderedStringVector(pointsTypeStr.c_str(), pointsStrings);
ASSERTL0(valid, "Unable to correctly parse the points types.");
for (vector<std::string>::size_type i = 0; i < pointsStrings.size(); i++)
{
valid = false;
for (unsigned int j = 0; j < LibUtilities::SIZE_PointsType; j++)
{
if (LibUtilities::kPointsTypeStr[j] == pointsStrings[i])
{
points.push_back((LibUtilities::PointsType) j);
valid = true;
break;
}
}
ASSERTL0(valid, std::string("Unable to correctly parse the points type: ").append(pointsStrings[i]).c_str());
}
const char *nPointsStr = expansion->Attribute("NUMPOINTS");
ASSERTL0(nPointsStr,"NUMPOINTS was not defined in EXPANSION section of input");
std::string numPointsStr = nPointsStr;
std::vector<unsigned int> numPoints;
valid = ParseUtils::GenerateOrderedVector(numPointsStr.c_str(), numPoints);
ASSERTL0(valid, "Unable to correctly parse the number of points.");
ASSERTL0(numPoints.size() == numPoints.size(),"information for num points does not match the number of basis");
for(int i = 0; i < basis.size(); ++i)
{
//Generate Basis key using information
const LibUtilities::PointsKey pkey(numPoints[i],points[i]);
basiskeyvec.push_back(LibUtilities::BasisKey(basis[i],numModes[i],pkey));
}
}
// Now have composite and basiskeys. Cycle through
// all composites for the geomShPtrs and set the modes
// and types for the elements contained in the element
// list.
CompositeMapIter compVecIter;
for (compVecIter = compositeVector.begin(); compVecIter != compositeVector.end(); ++compVecIter)
{
GeometryVectorIter geomVecIter;
for (geomVecIter = (compVecIter->second)->begin(); geomVecIter != (compVecIter->second)->end(); ++geomVecIter)
{
ExpansionMapIter x = expansionMap->find((*geomVecIter)->GetGlobalID());
ASSERTL0(x != expansionMap->end(), "Expansion not found!!");
if(useExpansionType)
{
(x->second)->m_basisKeyVector = MeshGraph::DefineBasisKeyFromExpansionType(*geomVecIter,expansion_type,num_modes);
}
else
{
ASSERTL0((*geomVecIter)->GetShapeDim() == basiskeyvec.size()," There is an incompatible expansion dimension with geometry dimension");
(x->second)->m_basisKeyVector = basiskeyvec;
}
}
}
expansion = expansion->NextSiblingElement("E");
}
}
else if(expType == "H")
{
int i;
ExpansionMapShPtr expansionMap;
while (expansion)
{
const char *fStr = expansion->Attribute("FIELDS");
std::vector<std::string> fieldStrings;
if(fStr) // extract other fields.
{
std::string fieldStr = fStr;
bool valid = ParseUtils::GenerateOrderedStringVector(fieldStr.c_str(),fieldStrings);
ASSERTL0(valid,"Unable to correctly parse the field string in ExpansionTypes.");
}
// check to see if m_expasionVectorShPtrMap has
// already been intiailised and if not intiailse
// vector.
if(m_expansionMapShPtrMap.count("DefaultVar") == 0) // no previous definitions
{
expansionMap = SetUpExpansionMap();
m_expansionMapShPtrMap["DefaultVar"] = expansionMap;
// make sure all fields in this search point
// to same expansion vector;
for(i = 0; i < fieldStrings.size(); ++i)
{
m_expansionMapShPtrMap[fieldStrings[i]] = expansionMap;
}
}
else // default variable is defined
{
if(fieldStrings.size()) // fields are defined
{
//see if field exists
if(m_expansionMapShPtrMap.count(fieldStrings[0]))
{
expansionMap = m_expansionMapShPtrMap.find(fieldStrings[0])->second;
}
else
{
expansionMap = SetUpExpansionMap();
// make sure all fields in this search point
// to same expansion vector;
for(i = 0; i < fieldStrings.size(); ++i)
{
if(m_expansionMapShPtrMap.count(fieldStrings[i]) == 0)
{
m_expansionMapShPtrMap[fieldStrings[i]] = expansionMap;
}
else
{
ASSERTL0(false,"Expansion vector for this field is already setup");
}
}
}
}
else // use default variable list
{
expansionMap = m_expansionMapShPtrMap.find("DefaultVar")->second;
}
}
/// Mandatory components...optional are to follow later.
std::string compositeStr = expansion->Attribute("COMPOSITE");
ASSERTL0(compositeStr.length() > 3, "COMPOSITE must be specified in expansion definition");
int beg = compositeStr.find_first_of("[");
int end = compositeStr.find_first_of("]");
std::string compositeListStr = compositeStr.substr(beg+1,end-beg-1);
CompositeMap compositeVector;
GetCompositeList(compositeListStr, compositeVector);
ExpansionType expansion_type_x = eNoExpansionType;
ExpansionType expansion_type_y = eNoExpansionType;
ExpansionType expansion_type_z = eNoExpansionType;
int num_modes_x = 0;
int num_modes_y = 0;
int num_modes_z = 0;
const char * tStr_x = expansion->Attribute("TYPE-X");
if(tStr_x) // use type string to define expansion
{
std::string typeStr = tStr_x;
const std::string* begStr = kExpansionTypeStr;
const std::string* endStr = kExpansionTypeStr+eExpansionTypeSize;
const std::string* expStr = std::find(begStr, endStr, typeStr);
ASSERTL0(expStr != endStr, "Invalid expansion type.");
expansion_type_x = (ExpansionType)(expStr - begStr);
const char *nStr = expansion->Attribute("NUMMODES-X");
ASSERTL0(nStr,"NUMMODES-X was not defined in EXPANSION section of input");
std::string nummodesStr = nStr;
// ASSERTL0(m_session,"Session should be defined to evaluate nummodes ");
if (m_session)
{
LibUtilities::Equation nummodesEqn(m_session, nummodesStr);
num_modes_x = (int) nummodesEqn.Evaluate();
}
else
{
num_modes_x = boost::lexical_cast<int>(nummodesStr);
}
}
const char * tStr_y = expansion->Attribute("TYPE-Y");
if(tStr_y) // use type string to define expansion
{
std::string typeStr = tStr_y;
const std::string* begStr = kExpansionTypeStr;
const std::string* endStr = kExpansionTypeStr+eExpansionTypeSize;
const std::string* expStr = std::find(begStr, endStr, typeStr);
ASSERTL0(expStr != endStr, "Invalid expansion type.");
expansion_type_y = (ExpansionType)(expStr - begStr);
const char *nStr = expansion->Attribute("NUMMODES-Y");
ASSERTL0(nStr,"NUMMODES-Y was not defined in EXPANSION section of input");
std::string nummodesStr = nStr;
// ASSERTL0(m_session,"Session should be defined to evaluate nummodes ");
if (m_session)
{
LibUtilities::Equation nummodesEqn(m_session, nummodesStr);
num_modes_y = (int) nummodesEqn.Evaluate();
}
else
{
num_modes_y = boost::lexical_cast<int>(nummodesStr);
}
}
const char * tStr_z = expansion->Attribute("TYPE-Z");
if(tStr_z) // use type string to define expansion
{
std::string typeStr = tStr_z;
const std::string* begStr = kExpansionTypeStr;
const std::string* endStr = kExpansionTypeStr+eExpansionTypeSize;
const std::string* expStr = std::find(begStr, endStr, typeStr);
ASSERTL0(expStr != endStr, "Invalid expansion type.");
expansion_type_z = (ExpansionType)(expStr - begStr);
const char *nStr = expansion->Attribute("NUMMODES-Z");
ASSERTL0(nStr,"NUMMODES-Z was not defined in EXPANSION section of input");
std::string nummodesStr = nStr;
// ASSERTL0(m_session,"Session should be defined to evaluate nummodes ");
if (m_session)
{
LibUtilities::Equation nummodesEqn(m_session, nummodesStr);
num_modes_z = (int) nummodesEqn.Evaluate();
}
else
{
num_modes_z = boost::lexical_cast<int>(nummodesStr);
}
}
CompositeMapIter compVecIter;
for (compVecIter = compositeVector.begin(); compVecIter != compositeVector.end(); ++compVecIter)
{
GeometryVectorIter geomVecIter;
for (geomVecIter = (compVecIter->second)->begin(); geomVecIter != (compVecIter->second)->end(); ++geomVecIter)
{
ExpansionMapIter expVecIter;
for (expVecIter = expansionMap->begin(); expVecIter != expansionMap->end(); ++expVecIter)
{
(expVecIter->second)->m_basisKeyVector = DefineBasisKeyFromExpansionTypeHomo(*geomVecIter,
expansion_type_x,
expansion_type_y,
expansion_type_z,
num_modes_x,
num_modes_y,
num_modes_z);
}
}
}
expansion = expansion->NextSiblingElement("H");
}
}
else if(expType == "ELEMENTS") // Reading a file with the expansion definition
{
std::vector<LibUtilities::FieldDefinitionsSharedPtr> fielddefs;
LibUtilities::FieldIO f(m_session->GetComm());
f.ImportFieldDefs(doc, fielddefs, true);
cout << " Number of elements: " << fielddefs.size() << endl;
SetExpansions(fielddefs);
}
else
{
ASSERTL0(false,"Expansion type not defined");
}
}
}
void Nektar::SpatialDomains::MeshGraph::ReadGeometry ( const std::string &  infilename)
virtual

Read will read the meshgraph vertices given a filename.

Reimplemented in Nektar::SpatialDomains::MeshGraph2D, Nektar::SpatialDomains::MeshGraph3D, and Nektar::SpatialDomains::MeshGraph1D.

Definition at line 229 of file MeshGraph.cpp.

References ASSERTL0.

Referenced by Read().

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

Read will read the meshgraph vertices given a TiXmlDocument.

Error value returned by TinyXML.

Reimplemented in Nektar::SpatialDomains::MeshGraph2D, Nektar::SpatialDomains::MeshGraph3D, and Nektar::SpatialDomains::MeshGraph1D.

Definition at line 247 of file MeshGraph.cpp.

References ASSERTL0, ASSERTL1, Nektar::LibUtilities::AnalyticExpressionEvaluator::DefineFunction(), Nektar::LibUtilities::AnalyticExpressionEvaluator::Evaluate(), m_meshDimension, m_meshPartitioned, m_partition, m_spaceDimension, and m_vertSet.

{
TiXmlHandle docHandle(&doc);
TiXmlElement* mesh = NULL;
TiXmlElement* master = NULL; // Master tag within which all data is contained.
int err; /// Error value returned by TinyXML.
master = doc.FirstChildElement("NEKTAR");
ASSERTL0(master, "Unable to find NEKTAR tag in file.");
// Find the Mesh tag and same the dim and space attributes
mesh = master->FirstChildElement("GEOMETRY");
ASSERTL0(mesh, "Unable to find GEOMETRY tag in file.");
TiXmlAttribute *attr = mesh->FirstAttribute();
// Initialize the mesh and space dimensions to 3 dimensions.
// We want to do this each time we read a file, so it should
// be done here and not just during class initialization.
while (attr)
{
std::string attrName(attr->Name());
if (attrName == "DIM")
{
err = attr->QueryIntValue(&m_meshDimension);
ASSERTL1(err==TIXML_SUCCESS, "Unable to read mesh dimension.");
}
else if (attrName == "SPACE")
{
err = attr->QueryIntValue(&m_spaceDimension);
ASSERTL1(err==TIXML_SUCCESS, "Unable to read space dimension.");
}
else if (attrName == "PARTITION")
{
err = attr->QueryIntValue(&m_partition);
ASSERTL1(err==TIXML_SUCCESS, "Unable to read partition.");
}
else
{
std::string errstr("Unknown attribute: ");
errstr += attrName;
ASSERTL1(false, errstr.c_str());
}
// Get the next attribute.
attr = attr->Next();
}
ASSERTL1(m_meshDimension<=m_spaceDimension, "Mesh dimension greater than space dimension");
// Now read the vertices
TiXmlElement* element = mesh->FirstChildElement("VERTEX");
ASSERTL0(element, "Unable to find mesh VERTEX tag in file.");
NekDouble xscale,yscale,zscale;
// check to see if any scaling parameters are in
// attributes and determine these values
LibUtilities::AnalyticExpressionEvaluator expEvaluator;
//LibUtilities::ExpressionEvaluator expEvaluator;
const char *xscal = element->Attribute("XSCALE");
if(!xscal)
{
xscale = 1.0;
}
else
{
std::string xscalstr = xscal;
int expr_id = expEvaluator.DefineFunction("",xscalstr);
xscale = expEvaluator.Evaluate(expr_id);
}
const char *yscal = element->Attribute("YSCALE");
if(!yscal)
{
yscale = 1.0;
}
else
{
std::string yscalstr = yscal;
int expr_id = expEvaluator.DefineFunction("",yscalstr);
yscale = expEvaluator.Evaluate(expr_id);
}
const char *zscal = element->Attribute("ZSCALE");
if(!zscal)
{
zscale = 1.0;
}
else
{
std::string zscalstr = zscal;
int expr_id = expEvaluator.DefineFunction("",zscalstr);
zscale = expEvaluator.Evaluate(expr_id);
}
NekDouble xmove,ymove,zmove;
// check to see if any moving parameters are in
// attributes and determine these values
//LibUtilities::ExpressionEvaluator expEvaluator;
const char *xmov = element->Attribute("XMOVE");
if(!xmov)
{
xmove = 0.0;
}
else
{
std::string xmovstr = xmov;
int expr_id = expEvaluator.DefineFunction("",xmovstr);
xmove = expEvaluator.Evaluate(expr_id);
}
const char *ymov = element->Attribute("YMOVE");
if(!ymov)
{
ymove = 0.0;
}
else
{
std::string ymovstr = ymov;
int expr_id = expEvaluator.DefineFunction("",ymovstr);
ymove = expEvaluator.Evaluate(expr_id);
}
const char *zmov = element->Attribute("ZMOVE");
if(!zmov)
{
zmove = 0.0;
}
else
{
std::string zmovstr = zmov;
int expr_id = expEvaluator.DefineFunction("",zmovstr);
zmove = expEvaluator.Evaluate(expr_id);
}
TiXmlElement *vertex = element->FirstChildElement("V");
int indx;
int nextVertexNumber = -1;
while (vertex)
{
nextVertexNumber++;
TiXmlAttribute *vertexAttr = vertex->FirstAttribute();
std::string attrName(vertexAttr->Name());
ASSERTL0(attrName == "ID", (std::string("Unknown attribute name: ") + attrName).c_str());
err = vertexAttr->QueryIntValue(&indx);
ASSERTL0(err == TIXML_SUCCESS, "Unable to read attribute ID.");
// Now read body of vertex
std::string vertexBodyStr;
TiXmlNode *vertexBody = vertex->FirstChild();
while (vertexBody)
{
// Accumulate all non-comment body data.
if (vertexBody->Type() == TiXmlNode::TEXT)
{
vertexBodyStr += vertexBody->ToText()->Value();
vertexBodyStr += " ";
}
vertexBody = vertexBody->NextSibling();
}
ASSERTL0(!vertexBodyStr.empty(), "Vertex definitions must contain vertex data.");
// Get vertex data from the data string.
NekDouble xval, yval, zval;
std::istringstream vertexDataStrm(vertexBodyStr.c_str());
try
{
while(!vertexDataStrm.fail())
{
vertexDataStrm >> xval >> yval >> zval;
xval = xval*xscale + xmove;
yval = yval*yscale + ymove;
zval = zval*zscale + zmove;
// Need to check it here because we may not be
// good after the read indicating that there
// was nothing to read.
if (!vertexDataStrm.fail())
{
PointGeomSharedPtr vert(MemoryManager<PointGeom>::AllocateSharedPtr(m_spaceDimension, indx, xval, yval, zval));
vert->SetGlobalID(indx);
m_vertSet[indx] = vert;
}
}
}
catch(...)
{
ASSERTL0(false, "Unable to read VERTEX data.");
}
vertex = vertex->NextSiblingElement("V");
}
}
void Nektar::SpatialDomains::MeshGraph::ReadGeometryInfo ( const std::string &  infilename)

Read geometric information from a file.

Read the geometry-related information from the given file. This information is located within the XML tree under <NEKTAR><GEOMETRY><GEOMINFO>.

Parameters
infilenameFilename of XML file.

Definition at line 469 of file MeshGraph.cpp.

References ASSERTL0.

{
TiXmlDocument doc(infilename);
bool loadOkay = doc.LoadFile();
std::stringstream errstr;
errstr << "Unable to load file: " << infilename << std::endl;
errstr << "Reason: " << doc.ErrorDesc() << std::endl;
errstr << "Position: Line " << doc.ErrorRow() << ", Column " << doc.ErrorCol() << std::endl;
ASSERTL0(loadOkay, errstr.str());
}
void Nektar::SpatialDomains::MeshGraph::ReadGeometryInfo ( TiXmlDocument &  doc)

Read geometric information from an XML document.

Read the geometry-related information from the given XML document. This information is located within the XML tree under <NEKTAR><GEOMETRY><GEOMINFO>.

Parameters
docXML document.

Definition at line 490 of file MeshGraph.cpp.

References ASSERTL0, Nektar::iterator, and m_geomInfo.

{
TiXmlElement *master = doc.FirstChildElement("NEKTAR");
ASSERTL0(master, "Unable to find NEKTAR tag in file.");
// Find the Expansions tag
TiXmlElement *geomTag = master->FirstChildElement("GEOMETRY");
ASSERTL0(geomTag, "Unable to find GEOMETRY tag in file.");
// See if we have GEOMINFO. If there is none, it's fine.
TiXmlElement *geomInfoTag = geomTag->FirstChildElement("GEOMINFO");
if (!geomInfoTag) return;
TiXmlElement *infoItem = geomInfoTag->FirstChildElement("I");
// Multiple nodes will only occur if there is a comment in between
// definitions.
while (infoItem)
{
std::string geomProperty = infoItem->Attribute("PROPERTY");
std::string geomValue = infoItem->Attribute("VALUE");
GeomInfoMap::iterator x = m_geomInfo.find(geomProperty);
ASSERTL0(x == m_geomInfo.end(),
"Property " + geomProperty + " already specified.");
m_geomInfo[geomProperty] = geomValue;
infoItem = infoItem->NextSiblingElement("I");
}
}
bool Nektar::SpatialDomains::MeshGraph::SameExpansions ( const std::string  var1,
const std::string  var2 
)
inline

Definition at line 506 of file MeshGraph.h.

References m_expansionMapShPtrMap.

{
ExpansionMapShPtr expVec1 = m_expansionMapShPtrMap.find(var1)->second;
ExpansionMapShPtr expVec2 = m_expansionMapShPtrMap.find(var2)->second;
if(expVec1.get() == expVec2.get())
{
return true;
}
return false;
}
void Nektar::SpatialDomains::MeshGraph::SetBasisKey ( LibUtilities::ShapeType  shape,
LibUtilities::BasisKeyVector keys,
std::string  var = "DefaultVar" 
)

Sets the basis key for all expansions of the given shape.

For each element of shape given by shape in field var, replace the current BasisKeyVector describing the expansion in each dimension, with the one provided by keys.

: Allow selection of elements through a CompositeVector, as well as by type.

Parameters
shapeThe shape of elements to be changed.
keysThe new basis vector to apply to those elements.

Definition at line 2565 of file MeshGraph.cpp.

References m_expansionMapShPtrMap.

{
ExpansionMapIter elemIter;
ExpansionMapShPtr expansionMap = m_expansionMapShPtrMap.find(var)->second;
for (elemIter = expansionMap->begin(); elemIter != expansionMap->end(); ++elemIter)
{
if ((elemIter->second)->m_geomShPtr->GetShapeType() == shape)
{
(elemIter->second)->m_basisKeyVector = keys;
}
}
}
void Nektar::SpatialDomains::MeshGraph::SetDomainRange ( NekDouble  xmin,
NekDouble  xmax,
NekDouble  ymin = NekConstants::kNekUnsetDouble,
NekDouble  ymax = NekConstants::kNekUnsetDouble,
NekDouble  zmin = NekConstants::kNekUnsetDouble,
NekDouble  zmax = NekConstants::kNekUnsetDouble 
)

Definition at line 1380 of file MeshGraph.cpp.

References Nektar::NekConstants::kNekUnsetDouble, m_domainRange, and Nektar::SpatialDomains::NullDomainRangeShPtr.

{
{
m_domainRange = MemoryManager<DomainRange>::AllocateSharedPtr();
m_domainRange->doXrange = true;
}
m_domainRange->xmin = xmin;
m_domainRange->xmax = xmax;
{
m_domainRange->doYrange = false;
}
else
{
m_domainRange->doYrange = true;
m_domainRange->ymin = ymin;
m_domainRange->ymax = ymax;
}
{
m_domainRange->doZrange = false;
}
else
{
m_domainRange->doZrange = true;
m_domainRange->zmin = zmin;
m_domainRange->zmax = zmax;
}
}
void Nektar::SpatialDomains::MeshGraph::SetExpansions ( std::vector< LibUtilities::FieldDefinitionsSharedPtr > &  fielddef)

Sets expansions given field definitions.

Definition at line 1751 of file MeshGraph.cpp.

References ASSERTL0, Nektar::LibUtilities::BasisTypeMap, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauMAlpha2Beta0, Nektar::LibUtilities::eHexahedron, Nektar::LibUtilities::ePrism, Nektar::LibUtilities::ePyramid, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eSegment, Nektar::LibUtilities::eTetrahedron, Nektar::LibUtilities::eTriangle, m_domain, m_expansionMapShPtrMap, m_hexGeoms, m_prismGeoms, m_pyrGeoms, m_quadGeoms, m_segGeoms, m_tetGeoms, and m_triGeoms.

Referenced by ReadExpansions().

{
int i, j, k, cnt, id;
ExpansionMapShPtr expansionMap;
// Loop over fields and determine unique fields string and
// declare whole expansion list
for(i = 0; i < fielddef.size(); ++i)
{
for(j = 0; j < fielddef[i]->m_fields.size(); ++j)
{
std::string field = fielddef[i]->m_fields[j];
if(m_expansionMapShPtrMap.count(field) == 0)
{
expansionMap = MemoryManager<ExpansionMap>::AllocateSharedPtr();
m_expansionMapShPtrMap[field] = expansionMap;
// check to see if DefaultVar also not set and if so assign it to this expansion
if(m_expansionMapShPtrMap.count("DefaultVar") == 0)
{
m_expansionMapShPtrMap["DefaultVar"] = expansionMap;
}
}
// loop over all elements in partition and set expansion
expansionMap = m_expansionMapShPtrMap.find(field)->second;
for(int d = 0; d < m_domain.size(); ++d)
{
CompositeMap::const_iterator compIter;
for (compIter = m_domain[d].begin();
compIter != m_domain[d].end(); ++compIter)
{
GeometryVector::const_iterator x;
for (x = compIter->second->begin();
x != compIter->second->end(); ++x)
{
ExpansionShPtr expansionElementShPtr =
MemoryManager<Expansion>::
AllocateSharedPtr(*x, def);
int id = (*x)->GetGlobalID();
(*expansionMap)[id] = expansionElementShPtr;
}
}
}
}
}
// loop over all elements find the geometry shared ptr and
// set up basiskey vector
for(i = 0; i < fielddef.size(); ++i)
{
cnt = 0;
std::vector<std::string> fields = fielddef[i]->m_fields;
std::vector<unsigned int> nmodes = fielddef[i]->m_numModes;
std::vector<LibUtilities::BasisType> basis = fielddef[i]->m_basis;
bool pointDef = fielddef[i]->m_pointsDef;
bool numPointDef = fielddef[i]->m_numPointsDef;
// Check points and numpoints
std::vector<unsigned int> npoints = fielddef[i]->m_numPoints;
std::vector<LibUtilities::PointsType> points = fielddef[i]->m_points;
bool UniOrder = fielddef[i]->m_uniOrder;
int check = 0;
for (j=0; j< basis.size(); ++j)
{
if ( (strcmp(LibUtilities::BasisTypeMap[basis[j]], "Modified_A") == 0) ||
(strcmp(LibUtilities::BasisTypeMap[basis[j]], "Modified_B") == 0) ||
(strcmp(LibUtilities::BasisTypeMap[basis[j]], "Modified_C") == 0) ||
(strcmp(LibUtilities::BasisTypeMap[basis[j]], "GLL_Lagrange") == 0) ||
(strcmp(LibUtilities::BasisTypeMap[basis[j]], "Gauss_Lagrange") == 0) ||
(strcmp(LibUtilities::BasisTypeMap[basis[j]], "Fourier") == 0) ||
(strcmp(LibUtilities::BasisTypeMap[basis[j]], "FourierSingleMode") == 0)||
(strcmp(LibUtilities::BasisTypeMap[basis[j]], "FourierHalfModeRe") == 0) ||
(strcmp(LibUtilities::BasisTypeMap[basis[j]], "FourierHalfModeIm") == 0))
{
check++;
}
}
if (check==basis.size())
{
for (j = 0; j < fielddef[i]->m_elementIDs.size(); ++j)
{
id = fielddef[i]->m_elementIDs[j];
switch (fielddef[i]->m_shapeType)
{
{
if(m_segGeoms.count(fielddef[i]->m_elementIDs[j]) == 0)
{
// skip element likely from parallel read
continue;
}
geom = m_segGeoms[fielddef[i]->m_elementIDs[j]];
LibUtilities::PointsKey pkey(nmodes[cnt]+1, LibUtilities::eGaussLobattoLegendre);
if(numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey1(npoints[cnt], points[0]);
pkey = pkey1;
}
else if(!numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey1(nmodes[cnt]+1, points[0]);
pkey = pkey1;
}
else if(numPointDef&&!pointDef)
{
const LibUtilities::PointsKey pkey1(npoints[cnt], LibUtilities::eGaussLobattoLegendre);
pkey = pkey1;
}
LibUtilities::BasisKey bkey(basis[0], nmodes[cnt], pkey);
if(!UniOrder)
{
cnt++;
}
bkeyvec.push_back(bkey);
}
break;
{
if(m_triGeoms.count(fielddef[i]->m_elementIDs[j]) == 0)
{
// skip element likely from parallel read
continue;
}
geom = m_triGeoms[fielddef[i]->m_elementIDs[j]];
LibUtilities::PointsKey pkey(nmodes[cnt]+1, LibUtilities::eGaussLobattoLegendre);
if(numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt], points[0]);
pkey = pkey2;
}
else if(!numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(nmodes[cnt]+1, points[0]);
pkey = pkey2;
}
else if(numPointDef&&!pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt], LibUtilities::eGaussLobattoLegendre);
pkey = pkey2;
}
LibUtilities::BasisKey bkey(basis[0], nmodes[cnt], pkey);
bkeyvec.push_back(bkey);
LibUtilities::PointsKey pkey1(nmodes[cnt+1], LibUtilities::eGaussRadauMAlpha1Beta0);
if(numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+1], points[1]);
pkey1 = pkey2;
}
else if(!numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(nmodes[cnt+1], points[1]);
pkey1 = pkey2;
}
else if(numPointDef&&!pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+1], LibUtilities::eGaussRadauMAlpha1Beta0);
pkey1 = pkey2;
}
LibUtilities::BasisKey bkey1(basis[1], nmodes[cnt+1], pkey1);
bkeyvec.push_back(bkey1);
if(!UniOrder)
{
cnt += 2;
}
}
break;
{
if(m_quadGeoms.count(fielddef[i]->m_elementIDs[j]) == 0)
{
// skip element likely from parallel read
continue;
}
geom = m_quadGeoms[fielddef[i]->m_elementIDs[j]];
for(int b = 0; b < 2; ++b)
{
LibUtilities::PointsKey pkey(nmodes[cnt+b]+1, LibUtilities::eGaussLobattoLegendre);
if(numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+b],points[b]);
pkey = pkey2;
}
else if(!numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(nmodes[cnt+b]+1,points[b]);
pkey = pkey2;
}
else if(numPointDef&&!pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+b],LibUtilities::eGaussLobattoLegendre);
pkey = pkey2;
}
LibUtilities::BasisKey bkey(basis[b],nmodes[cnt+b],pkey);
bkeyvec.push_back(bkey);
}
if(!UniOrder)
{
cnt += 2;
}
}
break;
{
k = fielddef[i]->m_elementIDs[j];
// allow for possibility that fielddef is
// larger than m_graph which can happen in
// parallel runs
if(m_tetGeoms.count(k) == 0)
{
continue;
}
geom = m_tetGeoms[k];
#if 0 //all gll
for(int b = 0; b < 3; ++b)
{
LibUtilities::PointsKey pkey(nmodes[cnt+b], LibUtilities::eGaussLobattoLegendre);
if(numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+b],points[b]);
pkey = pkey2;
}
else if(!numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(nmodes[cnt+b]+1,points[b]);
pkey = pkey2;
}
else if(numPointDef&&!pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+b],LibUtilities::eGaussLobattoLegendre);
pkey = pkey2;
}
LibUtilities::BasisKey bkey(basis[b],nmodes[cnt+b],pkey);
bkeyvec.push_back(bkey);
}
#else
{
LibUtilities::PointsKey pkey(nmodes[cnt], LibUtilities::eGaussLobattoLegendre);
if(numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt],points[0]);
pkey = pkey2;
}
else if(!numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(nmodes[cnt]+1,points[0]);
pkey = pkey2;
}
else if(numPointDef&&!pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt],LibUtilities::eGaussLobattoLegendre);
pkey = pkey2;
}
LibUtilities::BasisKey bkey(basis[0],nmodes[cnt],pkey);
bkeyvec.push_back(bkey);
}
{
LibUtilities::PointsKey pkey(nmodes[cnt+1], LibUtilities::eGaussRadauMAlpha1Beta0);
if(numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+1],points[1]);
pkey = pkey2;
}
else if(!numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(nmodes[cnt+1]+1,points[1]);
pkey = pkey2;
}
else if(numPointDef&&!pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+1],LibUtilities::eGaussRadauMAlpha1Beta0);
pkey = pkey2;
}
LibUtilities::BasisKey bkey(basis[1],nmodes[cnt+1],pkey);
bkeyvec.push_back(bkey);
}
{
LibUtilities::PointsKey pkey(nmodes[cnt+2], LibUtilities::eGaussRadauMAlpha2Beta0);
if(numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+2],points[2]);
pkey = pkey2;
}
else if(!numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(nmodes[cnt+2]+1,points[2]);
pkey = pkey2;
}
else if(numPointDef&&!pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+2],LibUtilities::eGaussRadauMAlpha1Beta0);
pkey = pkey2;
}
LibUtilities::BasisKey bkey(basis[2],nmodes[cnt+2],pkey);
bkeyvec.push_back(bkey);
}
#endif
if(!UniOrder)
{
cnt += 3;
}
}
break;
{
k = fielddef[i]->m_elementIDs[j];
if(m_prismGeoms.count(k) == 0)
{
continue;
}
geom = m_prismGeoms[k];
#if 0 // all GLL
for(int b = 0; b < 3; ++b)
{
LibUtilities::PointsKey pkey(nmodes[cnt+b],LibUtilities::eGaussLobattoLegendre);
if(numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+b],points[b]);
pkey = pkey2;
}
else if(!numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(nmodes[cnt+b]+1,points[b]);
pkey = pkey2;
}
else if(numPointDef&&!pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+b],LibUtilities::eGaussLobattoLegendre);
pkey = pkey2;
}
LibUtilities::BasisKey bkey(basis[b],nmodes[cnt+b],pkey);
bkeyvec.push_back(bkey);
}
#else
for(int b = 0; b < 2; ++b)
{
LibUtilities::PointsKey pkey(nmodes[cnt+b],LibUtilities::eGaussLobattoLegendre);
if(numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+b],points[b]);
pkey = pkey2;
}
else if(!numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(nmodes[cnt+b]+1,points[b]);
pkey = pkey2;
}
else if(numPointDef&&!pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+b],LibUtilities::eGaussLobattoLegendre);
pkey = pkey2;
}
LibUtilities::BasisKey bkey(basis[b],nmodes[cnt+b],pkey);
bkeyvec.push_back(bkey);
}
{
LibUtilities::PointsKey pkey(nmodes[cnt+2],LibUtilities::eGaussRadauMAlpha1Beta0);
if(numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+2],points[2]);
pkey = pkey2;
}
else if(!numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(nmodes[cnt+2]+1,points[2]);
pkey = pkey2;
}
else if(numPointDef&&!pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+2],LibUtilities::eGaussLobattoLegendre);
pkey = pkey2;
}
LibUtilities::BasisKey bkey(basis[2],nmodes[cnt+2],pkey);
bkeyvec.push_back(bkey);
}
#endif
if(!UniOrder)
{
cnt += 3;
}
}
break;
{
k = fielddef[i]->m_elementIDs[j];
ASSERTL0(m_pyrGeoms.find(k) != m_pyrGeoms.end(),
"Failed to find geometry with same global id");
geom = m_pyrGeoms[k];
for(int b = 0; b < 3; ++b)
{
LibUtilities::PointsKey pkey(nmodes[cnt+b],points[b]);
if(numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+b],points[b]);
pkey = pkey2;
}
else if(!numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(nmodes[cnt+b]+1,points[b]);
pkey = pkey2;
}
else if(numPointDef&&!pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+b],LibUtilities::eGaussLobattoLegendre);
pkey = pkey2;
}
LibUtilities::BasisKey bkey(basis[b],nmodes[cnt+b],pkey);
bkeyvec.push_back(bkey);
}
if(!UniOrder)
{
cnt += 3;
}
}
break;
{
k = fielddef[i]->m_elementIDs[j];
if(m_hexGeoms.count(k) == 0)
{
continue;
}
geom = m_hexGeoms[k];
for(int b = 0; b < 3; ++b)
{
LibUtilities::PointsKey pkey(nmodes[cnt+b],LibUtilities::eGaussLobattoLegendre);
if(numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+b],points[b]);
pkey = pkey2;
}
else if(!numPointDef&&pointDef)
{
const LibUtilities::PointsKey pkey2(nmodes[cnt+b]+1,points[b]);
pkey = pkey2;
}
else if(numPointDef&&!pointDef)
{
const LibUtilities::PointsKey pkey2(npoints[cnt+b],LibUtilities::eGaussLobattoLegendre);
pkey = pkey2;
}
LibUtilities::BasisKey bkey(basis[b],nmodes[cnt+b],pkey);
bkeyvec.push_back(bkey);
}
if(!UniOrder)
{
cnt += 3;
}
}
break;
default:
ASSERTL0(false,"Need to set up for pyramid and prism 3D Expansions");
break;
}
for(k = 0; k < fields.size(); ++k)
{
expansionMap = m_expansionMapShPtrMap.find(fields[k])->second;
if((*expansionMap).find(id) != (*expansionMap).end())
{
(*expansionMap)[id]->m_geomShPtr = geom;
(*expansionMap)[id]->m_basisKeyVector = bkeyvec;
}
}
}
}
else
{
ASSERTL0(false,"Need to set up for non Modified basis");
}
}
}
void Nektar::SpatialDomains::MeshGraph::SetExpansions ( std::vector< LibUtilities::FieldDefinitionsSharedPtr > &  fielddef,
std::vector< std::vector< LibUtilities::PointsType > > &  pointstype 
)

Sets expansions given field definition, quadrature points.

Definition at line 2286 of file MeshGraph.cpp.

References ASSERTL0, Nektar::LibUtilities::eHexahedron, Nektar::LibUtilities::ePrism, Nektar::LibUtilities::ePyramid, Nektar::LibUtilities::eQuadrilateral, Nektar::LibUtilities::eSegment, Nektar::LibUtilities::eTetrahedron, Nektar::LibUtilities::eTriangle, m_expansionMapShPtrMap, m_hexGeoms, m_prismGeoms, m_pyrGeoms, m_quadGeoms, m_segGeoms, m_tetGeoms, and m_triGeoms.

{
int i,j,k,g,h,cnt,id;
ExpansionMapShPtr expansionMap;
// Loop over fields and determine unique fields string and
// declare whole expansion list
for(i = 0; i < fielddef.size(); ++i)
{
for(j = 0; j < fielddef[i]->m_fields.size(); ++j)
{
std::string field = fielddef[i]->m_fields[j];
if(m_expansionMapShPtrMap.count(field) == 0)
{
expansionMap = MemoryManager<ExpansionMap>::AllocateSharedPtr();
m_expansionMapShPtrMap[field] = expansionMap;
// check to see if DefaultVar also not set and if so assign it to this expansion
if(m_expansionMapShPtrMap.count("DefaultVar") == 0)
{
m_expansionMapShPtrMap["DefaultVar"] = expansionMap;
}
// loop over all elements and set expansion
for(k = 0; k < fielddef.size(); ++k)
{
for(h = 0; h < fielddef[k]->m_fields.size(); ++h)
{
if(fielddef[k]->m_fields[h] == field)
{
expansionMap = m_expansionMapShPtrMap.find(field)->second;
for(g = 0; g < fielddef[k]->m_elementIDs.size(); ++g)
{
ExpansionShPtr tmpexp =
MemoryManager<Expansion>::AllocateSharedPtr(geom, def);
(*expansionMap)[fielddef[k]->m_elementIDs[g]] = tmpexp;
}
}
}
}
}
}
}
// loop over all elements find the geometry shared ptr and
// set up basiskey vector
for(i = 0; i < fielddef.size(); ++i)
{
cnt = 0;
std::vector<std::string> fields = fielddef[i]->m_fields;
std::vector<unsigned int> nmodes = fielddef[i]->m_numModes;
std::vector<LibUtilities::BasisType> basis = fielddef[i]->m_basis;
bool UniOrder = fielddef[i]->m_uniOrder;
for(j = 0; j < fielddef[i]->m_elementIDs.size(); ++j)
{
id = fielddef[i]->m_elementIDs[j];
switch(fielddef[i]->m_shapeType)
{
{
k = fielddef[i]->m_elementIDs[j];
ASSERTL0(m_segGeoms.find(k) != m_segGeoms.end(),
"Failed to find geometry with same global id.");
geom = m_segGeoms[k];
const LibUtilities::PointsKey pkey(nmodes[cnt], pointstype[i][0]);
LibUtilities::BasisKey bkey(basis[0], nmodes[cnt], pkey);
if(!UniOrder)
{
cnt++;
}
bkeyvec.push_back(bkey);
}
break;
{
k = fielddef[i]->m_elementIDs[j];
ASSERTL0(m_triGeoms.find(k) != m_triGeoms.end(),
"Failed to find geometry with same global id.");
geom = m_triGeoms[k];
for(int b = 0; b < 2; ++b)
{
const LibUtilities::PointsKey pkey(nmodes[cnt+b],pointstype[i][b]);
LibUtilities::BasisKey bkey(basis[b],nmodes[cnt+b],pkey);
bkeyvec.push_back(bkey);
}
if(!UniOrder)
{
cnt += 2;
}
}
break;
{
k = fielddef[i]->m_elementIDs[j];
ASSERTL0(m_quadGeoms.find(k) != m_quadGeoms.end(),
"Failed to find geometry with same global id");
geom = m_quadGeoms[k];
for(int b = 0; b < 2; ++b)
{
const LibUtilities::PointsKey pkey(nmodes[cnt+b],pointstype[i][b]);
LibUtilities::BasisKey bkey(basis[b],nmodes[cnt+b],pkey);
bkeyvec.push_back(bkey);
}
if(!UniOrder)
{
cnt += 2;
}
}
break;
{
k = fielddef[i]->m_elementIDs[j];
ASSERTL0(m_tetGeoms.find(k) != m_tetGeoms.end(),
"Failed to find geometry with same global id");
geom = m_tetGeoms[k];
for(int b = 0; b < 3; ++b)
{
const LibUtilities::PointsKey pkey(nmodes[cnt+b],pointstype[i][b]);
LibUtilities::BasisKey bkey(basis[b],nmodes[cnt+b],pkey);
bkeyvec.push_back(bkey);
}
if(!UniOrder)
{
cnt += 2;
}
}
break;
{
k = fielddef[i]->m_elementIDs[j];
ASSERTL0(m_pyrGeoms.find(k) != m_pyrGeoms.end(),
"Failed to find geometry with same global id");
geom = m_pyrGeoms[k];
for(int b = 0; b < 3; ++b)
{
const LibUtilities::PointsKey pkey(nmodes[cnt+b],pointstype[i][b]);
LibUtilities::BasisKey bkey(basis[b],nmodes[cnt+b],pkey);
bkeyvec.push_back(bkey);
}
if(!UniOrder)
{
cnt += 2;
}
}
break;
{
k = fielddef[i]->m_elementIDs[j];
"Failed to find geometry with same global id");
geom = m_prismGeoms[k];
for(int b = 0; b < 3; ++b)
{
const LibUtilities::PointsKey pkey(nmodes[cnt+b],pointstype[i][b]);
LibUtilities::BasisKey bkey(basis[b],nmodes[cnt+b],pkey);
bkeyvec.push_back(bkey);
}
if(!UniOrder)
{
cnt += 2;
}
}
break;
{
k = fielddef[i]->m_elementIDs[j];
ASSERTL0(m_hexGeoms.find(k) != m_hexGeoms.end(),
"Failed to find geometry with same global id");
geom = m_hexGeoms[k];
for(int b = 0; b < 3; ++b)
{
const LibUtilities::PointsKey pkey(nmodes[cnt+b],pointstype[i][b]);
LibUtilities::BasisKey bkey(basis[b],nmodes[cnt+b],pkey);
bkeyvec.push_back(bkey);
}
if(!UniOrder)
{
cnt += 2;
}
}
break;
default:
ASSERTL0(false,"Need to set up for pyramid and prism 3D Expansions");
break;
}
for(k = 0; k < fields.size(); ++k)
{
expansionMap = m_expansionMapShPtrMap.find(fields[k])->second;
if((*expansionMap).find(id) != (*expansionMap).end())
{
(*expansionMap)[id]->m_geomShPtr = geom;
(*expansionMap)[id]->m_basisKeyVector = bkeyvec;
}
}
}
}
}
void Nektar::SpatialDomains::MeshGraph::SetExpansions ( const std::string  variable,
ExpansionMapShPtr exp 
)
inline

This function sets the expansion #exp in map with entry #variable.

Definition at line 490 of file MeshGraph.h.

References ASSERTL0, and m_expansionMapShPtrMap.

{
if(m_expansionMapShPtrMap.count(variable) != 0)
{
ASSERTL0(false,(std::string("Expansion field is already set for variable ") + variable).c_str());
}
else
{
m_expansionMapShPtrMap[variable] = exp;
}
}
void Nektar::SpatialDomains::MeshGraph::SetExpansionsToEvenlySpacedPoints ( int  npoints = 0)

Sets expansions to have equispaced points.

Reset all points keys to have equispaced points with optional arguemtn of npoints which redefines how many points are to be used.

Definition at line 2513 of file MeshGraph.cpp.

References Nektar::LibUtilities::ePolyEvenlySpaced, Nektar::LibUtilities::BasisKey::GetBasisType(), Nektar::LibUtilities::BasisKey::GetNumModes(), m_expansionMapShPtrMap, and npts.

{
// iterate over all defined expansions
for(it = m_expansionMapShPtrMap.begin(); it != m_expansionMapShPtrMap.end(); ++it)
{
for(expIt = it->second->begin(); expIt != it->second->end(); ++expIt)
{
for(int i = 0; i < expIt->second->m_basisKeyVector.size(); ++i)
{
LibUtilities::BasisKey bkeyold = expIt->second->m_basisKeyVector[i];
int npts;
if(npoints) // use input
{
npts = npoints;
}
else
{
npts = bkeyold.GetNumModes();
}
const LibUtilities::PointsKey pkey(npts,LibUtilities::ePolyEvenlySpaced);
LibUtilities::BasisKey bkeynew(bkeyold.GetBasisType(),bkeyold.GetNumModes(), pkey);
expIt->second->m_basisKeyVector[i] = bkeynew;
}
}
}
}
ExpansionMapShPtr Nektar::SpatialDomains::MeshGraph::SetUpExpansionMap ( void  )
protected

Generate a single vector of Expansion structs mapping global element ID to a corresponding Geometry shared pointer and basis key.

Expansion map ensures elements which appear in multiple composites within the domain are only listed once.

Definition at line 3503 of file MeshGraph.cpp.

References m_domain.

Referenced by ReadExpansions().

{
ExpansionMapShPtr returnval;
returnval = MemoryManager<ExpansionMap>::AllocateSharedPtr();
for(int d = 0; d < m_domain.size(); ++d)
{
CompositeMap::const_iterator compIter;
for (compIter = m_domain[d].begin(); compIter != m_domain[d].end(); ++compIter)
{
GeometryVector::const_iterator x;
for (x = compIter->second->begin(); x != compIter->second->end(); ++x)
{
ExpansionShPtr expansionElementShPtr =
MemoryManager<Expansion>::AllocateSharedPtr(*x, def);
int id = (*x)->GetGlobalID();
(*returnval)[id] = expansionElementShPtr;
}
}
}
return returnval;
}

Member Data Documentation

CurveVector Nektar::SpatialDomains::MeshGraph::m_curvedEdges
protected

Definition at line 389 of file MeshGraph.h.

Referenced by GetCurvedEdges(), ReadCurves(), Nektar::SpatialDomains::MeshGraph3D::ReadEdges(), and Nektar::SpatialDomains::MeshGraph2D::ReadEdges().

CurveVector Nektar::SpatialDomains::MeshGraph::m_curvedFaces
protected

Definition at line 390 of file MeshGraph.h.

Referenced by GetCurvedFaces(), ReadCurves(), Nektar::SpatialDomains::MeshGraph2D::ReadElements(), and Nektar::SpatialDomains::MeshGraph3D::ReadFaces().

std::vector<CompositeMap> Nektar::SpatialDomains::MeshGraph::m_domain
protected

Definition at line 407 of file MeshGraph.h.

Referenced by GetDomain(), Nektar::SpatialDomains::MeshGraph2D::GetElementsFromEdge(), ReadDomain(), SetExpansions(), and SetUpExpansionMap().

DomainRangeShPtr Nektar::SpatialDomains::MeshGraph::m_domainRange
protected

Definition at line 408 of file MeshGraph.h.

Referenced by CheckRange(), and SetDomainRange().

ExpansionMapShPtrMap Nektar::SpatialDomains::MeshGraph::m_expansionMapShPtrMap
protected

Definition at line 410 of file MeshGraph.h.

Referenced by GetExpansion(), GetExpansions(), ReadExpansions(), SameExpansions(), SetBasisKey(), SetExpansions(), and SetExpansionsToEvenlySpacedPoints().

GeomInfoMap Nektar::SpatialDomains::MeshGraph::m_geomInfo
protected

Definition at line 412 of file MeshGraph.h.

Referenced by CheckForGeomInfo(), GetGeomInfo(), and ReadGeometryInfo().

HexGeomMap Nektar::SpatialDomains::MeshGraph::m_hexGeoms
protected

Definition at line 399 of file MeshGraph.h.

Referenced by AddHexahedron(), GetAllHexGeoms(), Nektar::SpatialDomains::MeshGraph3D::ReadElements(), Nektar::SpatialDomains::MeshGraph3D::ResolveGeomRef(), and SetExpansions().

InterfaceCompList Nektar::SpatialDomains::MeshGraph::m_iComps
protected

Definition at line 387 of file MeshGraph.h.

CompositeMap Nektar::SpatialDomains::MeshGraph::m_meshComposites
protected

Definition at line 406 of file MeshGraph.h.

Referenced by GetComposite(), GetCompositeItem(), GetCompositeList(), GetComposites(), Nektar::SpatialDomains::MeshGraph3D::GetNumCompositeItems(), Nektar::SpatialDomains::MeshGraph2D::GetNumCompositeItems(), Nektar::SpatialDomains::MeshGraph3D::GetNumComposites(), Nektar::SpatialDomains::MeshGraph2D::GetNumComposites(), Nektar::SpatialDomains::MeshGraph1D::ReadComposites(), Nektar::SpatialDomains::MeshGraph3D::ReadComposites(), and Nektar::SpatialDomains::MeshGraph2D::ReadComposites().

int Nektar::SpatialDomains::MeshGraph::m_meshDimension
protected

Definition at line 401 of file MeshGraph.h.

Referenced by GetMeshDimension(), ReadCurves(), and ReadGeometry().

bool Nektar::SpatialDomains::MeshGraph::m_meshPartitioned
protected

Definition at line 404 of file MeshGraph.h.

Referenced by GetCompositeList(), and ReadGeometry().

int Nektar::SpatialDomains::MeshGraph::m_partition
protected

Definition at line 403 of file MeshGraph.h.

Referenced by ReadGeometry().

PrismGeomMap Nektar::SpatialDomains::MeshGraph::m_prismGeoms
protected

Definition at line 398 of file MeshGraph.h.

Referenced by AddPrism(), GetAllPrismGeoms(), Nektar::SpatialDomains::MeshGraph3D::ReadElements(), Nektar::SpatialDomains::MeshGraph3D::ResolveGeomRef(), and SetExpansions().

PyrGeomMap Nektar::SpatialDomains::MeshGraph::m_pyrGeoms
protected

Definition at line 397 of file MeshGraph.h.

Referenced by AddPyramid(), GetAllPyrGeoms(), Nektar::SpatialDomains::MeshGraph3D::ReadElements(), Nektar::SpatialDomains::MeshGraph3D::ResolveGeomRef(), and SetExpansions().

QuadGeomMap Nektar::SpatialDomains::MeshGraph::m_quadGeoms
protected

Definition at line 395 of file MeshGraph.h.

Referenced by AddQuadrilateral(), GetAllQuadGeoms(), Nektar::SpatialDomains::MeshGraph3D::GetCartesianEorientFromElmt(), Nektar::SpatialDomains::MeshGraph2D::GetCartesianEorientFromElmt(), Nektar::SpatialDomains::MeshGraph3D::GetEidFromElmt(), Nektar::SpatialDomains::MeshGraph2D::GetEidFromElmt(), Nektar::SpatialDomains::MeshGraph3D::GetEorientFromElmt(), Nektar::SpatialDomains::MeshGraph2D::GetEorientFromElmt(), Nektar::SpatialDomains::MeshGraph3D::GetGeometry2D(), Nektar::SpatialDomains::MeshGraph2D::GetQuadgeoms(), Nektar::SpatialDomains::MeshGraph3D::GetQuadgeoms(), Nektar::SpatialDomains::MeshGraph3D::GetVidFromElmt(), Nektar::SpatialDomains::MeshGraph2D::GetVidFromElmt(), Nektar::SpatialDomains::MeshGraph2D::ReadElements(), Nektar::SpatialDomains::MeshGraph3D::ReadFaces(), Nektar::SpatialDomains::MeshGraph3D::ResolveGeomRef(), Nektar::SpatialDomains::MeshGraph2D::ResolveGeomRef(), and SetExpansions().

SegGeomMap Nektar::SpatialDomains::MeshGraph::m_segGeoms
protected

Definition at line 392 of file MeshGraph.h.

Referenced by AddEdge(), GetAllSegGeoms(), GetEdge(), Nektar::SpatialDomains::MeshGraph3D::GetNseggeoms(), Nektar::SpatialDomains::MeshGraph2D::GetNseggeoms(), Nektar::SpatialDomains::MeshGraph3D::GetSegGeom(), Nektar::SpatialDomains::MeshGraph2D::GetSegGeom(), Nektar::SpatialDomains::MeshGraph1D::GetSeggeoms(), Nektar::SpatialDomains::MeshGraph1D::GetVidFromElmt(), Nektar::SpatialDomains::MeshGraph3D::ReadEdges(), Nektar::SpatialDomains::MeshGraph2D::ReadEdges(), Nektar::SpatialDomains::MeshGraph1D::ReadElements(), Nektar::SpatialDomains::MeshGraph1D::ResolveGeomRef(), Nektar::SpatialDomains::MeshGraph3D::ResolveGeomRef(), Nektar::SpatialDomains::MeshGraph2D::ResolveGeomRef(), and SetExpansions().

LibUtilities::SessionReaderSharedPtr Nektar::SpatialDomains::MeshGraph::m_session
protected

Definition at line 385 of file MeshGraph.h.

Referenced by ReadExpansions().

int Nektar::SpatialDomains::MeshGraph::m_spaceDimension
protected

Definition at line 402 of file MeshGraph.h.

Referenced by AddEdge(), AddVertex(), GetSpaceDimension(), Nektar::SpatialDomains::MeshGraph3D::ReadEdges(), Nektar::SpatialDomains::MeshGraph2D::ReadEdges(), and ReadGeometry().

TetGeomMap Nektar::SpatialDomains::MeshGraph::m_tetGeoms
protected

Definition at line 396 of file MeshGraph.h.

Referenced by AddTetrahedron(), GetAllTetGeoms(), Nektar::SpatialDomains::MeshGraph3D::ReadElements(), Nektar::SpatialDomains::MeshGraph3D::ResolveGeomRef(), and SetExpansions().

TriGeomMap Nektar::SpatialDomains::MeshGraph::m_triGeoms
protected

Definition at line 394 of file MeshGraph.h.

Referenced by AddTriangle(), GetAllTriGeoms(), Nektar::SpatialDomains::MeshGraph3D::GetCartesianEorientFromElmt(), Nektar::SpatialDomains::MeshGraph2D::GetCartesianEorientFromElmt(), Nektar::SpatialDomains::MeshGraph3D::GetEidFromElmt(), Nektar::SpatialDomains::MeshGraph2D::GetEidFromElmt(), Nektar::SpatialDomains::MeshGraph3D::GetEorientFromElmt(), Nektar::SpatialDomains::MeshGraph2D::GetEorientFromElmt(), Nektar::SpatialDomains::MeshGraph3D::GetGeometry2D(), Nektar::SpatialDomains::MeshGraph2D::GetTrigeoms(), Nektar::SpatialDomains::MeshGraph3D::GetTrigeoms(), Nektar::SpatialDomains::MeshGraph3D::GetVidFromElmt(), Nektar::SpatialDomains::MeshGraph2D::GetVidFromElmt(), Nektar::SpatialDomains::MeshGraph2D::ReadElements(), Nektar::SpatialDomains::MeshGraph3D::ReadFaces(), Nektar::SpatialDomains::MeshGraph3D::ResolveGeomRef(), Nektar::SpatialDomains::MeshGraph2D::ResolveGeomRef(), and SetExpansions().

PointGeomMap Nektar::SpatialDomains::MeshGraph::m_vertSet
protected

Definition at line 386 of file MeshGraph.h.

Referenced by AddVertex(), GetNvertices(), GetVertex(), ReadGeometry(), Nektar::SpatialDomains::MeshGraph1D::ResolveGeomRef(), Nektar::SpatialDomains::MeshGraph3D::ResolveGeomRef(), and Nektar::SpatialDomains::MeshGraph2D::ResolveGeomRef().