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

A 3-dimensional five-faced element (2 triangles, 3 quadrilaterals). More...

#include <MeshElements.h>

Inheritance diagram for Nektar::Utilities::Prism:
Inheritance graph
[legend]
Collaboration diagram for Nektar::Utilities::Prism:
Collaboration graph
[legend]

Public Member Functions

 Prism (ElmtConfig pConf, std::vector< NodeSharedPtr > pNodeList, std::vector< int > pTagList)
 Create a prism element.
 Prism (const Prism &pSrc)
virtual ~Prism ()
virtual
SpatialDomains::GeometrySharedPtr 
GetGeom (int coordDim)
 Generate a Nektar++ geometry object for this element.
virtual void Complete (int order)
 
- Public Member Functions inherited from Nektar::Utilities::Element
 Element (ElmtConfig pConf, unsigned int pNumNodes, unsigned int pGotNodes)
unsigned int GetId () const
 Returns the ID of the element (or associated edge or face for boundary elements).
unsigned int GetDim () const
 Returns the expansion dimension of the element.
ElmtConfig GetConf () const
 Returns the configuration of the element.
std::string GetTag () const
 Returns the tag which defines the element shape.
NodeSharedPtr GetVertex (unsigned int i) const
 Access a vertex node.
EdgeSharedPtr GetEdge (unsigned int i) const
 Access an edge.
FaceSharedPtr GetFace (unsigned int i) const
 Access a face.
std::vector< NodeSharedPtrGetVertexList () const
 Access the list of vertex nodes.
std::vector< EdgeSharedPtrGetEdgeList () const
 Access the list of edges.
std::vector< FaceSharedPtrGetFaceList () const
 Access the list of faces.
std::vector< NodeSharedPtrGetVolumeNodes () const
 Access the list of volume nodes.
void SetVolumeNodes (std::vector< NodeSharedPtr > &nodes)
LibUtilities::PointsType GetCurveType () const
void SetCurveType (LibUtilities::PointsType cT)
unsigned int GetNodeCount () const
 Returns the total number of nodes (vertices, edge nodes and face nodes and volume nodes).
std::vector< int > GetTagList () const
 Access the list of tags associated with this element.
unsigned int GetVertexCount () const
 Returns the number of vertices.
unsigned int GetEdgeCount () const
 Returns the number of edges.
unsigned int GetFaceCount () const
 Returns the number of faces.
void SetId (unsigned int p)
 Change the ID of the element.
void SetVertex (unsigned int p, NodeSharedPtr pNew)
 Replace a vertex with another vertex object.
void SetEdge (unsigned int p, EdgeSharedPtr pNew)
 Replace an edge with another edge object.
void SetFace (unsigned int p, FaceSharedPtr pNew)
 Replace a face with another face object.
void SetEdgeLink (EdgeSharedPtr pLink)
 Set a correspondence between this element and an edge (2D boundary element).
EdgeSharedPtr GetEdgeLink ()
 Get correspondence between this element and an edge.
void SetFaceLink (FaceSharedPtr pLink)
 Set a correspondence between this element and a face (3D boundary element).
FaceSharedPtr GetFaceLink ()
 Get correspondence between this element and a face.
void SetBoundaryLink (int i, int j)
 Set a correspondence between edge or face i and its representative boundary element m->element[expDim-1][j].
int GetBoundaryLink (int i)
 Get the location of the boundary face/edge i for this element.
void SetTagList (const std::vector< int > &tags)
 Set the list of tags associated with this element.
virtual std::string GetXmlString () const
 Generate a list of vertices (1D), edges (2D), or faces (3D).
std::vector< NodeSharedPtrtensorNodeOrdering (const std::vector< NodeSharedPtr > &nodes, int n) const
 Reorders Gmsh ordered nodes into a row-by-row ordering required for Nektar++ curve tags.
std::string GetXmlCurveString () const
 Generates a string listing the coordinates of all nodes associated with this element.
int GetMaxOrder ()
 Obtain the order of an element by looking at edges.
void Print ()

Static Public Member Functions

static ElementSharedPtr create (ElmtConfig pConf, std::vector< NodeSharedPtr > pNodeList, std::vector< int > pTagList)
 Creates an instance of this class.
static unsigned int GetNumNodes (ElmtConfig pConf)
 Return the number of nodes defining a prism.

Public Attributes

unsigned int m_orientation

Static Public Attributes

static LibUtilities::ShapeType m_type
 Element type.

Protected Member Functions

void OrientPrism ()
 Orient prism to align degenerate vertices.

Additional Inherited Members

- Protected Attributes inherited from Nektar::Utilities::Element
unsigned int m_id
 ID of the element.
unsigned int m_dim
 Dimension of the element.
ElmtConfig m_conf
 Contains configuration of the element.
std::string m_tag
 Tag character describing the element.
std::vector< int > m_taglist
 List of integers specifying properties of the element.
std::vector< NodeSharedPtrm_vertex
 List of element vertex nodes.
std::vector< EdgeSharedPtrm_edge
 List of element edges.
std::vector< FaceSharedPtrm_face
 List of element faces.
std::vector< NodeSharedPtrm_volumeNodes
 List of element volume nodes.
LibUtilities::PointsType m_curveType
 Volume curve type.
EdgeSharedPtr m_edgeLink
 Pointer to the corresponding edge if element is a 2D boundary.
FaceSharedPtr m_faceLink
 Pointer to the corresponding face if element is a 3D boundary.
std::map< int, int > m_boundaryLinks
 Array mapping faces/edges to the location of the appropriate boundary elements in m->element.
SpatialDomains::GeometrySharedPtr m_geom
 Nektar++ geometry object for this element.

Detailed Description

A 3-dimensional five-faced element (2 triangles, 3 quadrilaterals).

Definition at line 1422 of file MeshElements.h.

Constructor & Destructor Documentation

Nektar::Utilities::Prism::Prism ( ElmtConfig  pConf,
std::vector< NodeSharedPtr pNodeList,
std::vector< int >  pTagList 
)

Create a prism element.

Definition at line 1350 of file MeshElements.cpp.

References Nektar::iterator, Nektar::Utilities::Element::m_conf, Nektar::Utilities::Element::m_dim, Nektar::Utilities::Element::m_edge, Nektar::Utilities::ElmtConfig::m_edgeCurveType, Nektar::Utilities::Element::m_face, Nektar::Utilities::ElmtConfig::m_faceCurveType, Nektar::Utilities::ElmtConfig::m_faceNodes, Nektar::Utilities::ElmtConfig::m_order, Nektar::Utilities::ElmtConfig::m_reorient, Nektar::Utilities::Element::m_tag, Nektar::Utilities::Element::m_taglist, Nektar::Utilities::Element::m_vertex, and OrientPrism().

: Element(pConf, GetNumNodes(pConf), pNodeList.size())
{
m_tag = "R";
m_dim = 3;
m_taglist = pTagList;
int n = m_conf.m_order-1;
// Create a map to relate edge nodes to a pair of vertices
// defining an edge. This is based on the ordering produced by
// gmsh.
map<pair<int,int>, int> edgeNodeMap;
map<pair<int,int>, int>::iterator it;
// This edge-node map is based on Nektar++ ordering.
edgeNodeMap[pair<int,int>(1,2)] = 7;
edgeNodeMap[pair<int,int>(2,3)] = 7 + n;
edgeNodeMap[pair<int,int>(4,3)] = 7 + 2*n;
edgeNodeMap[pair<int,int>(1,4)] = 7 + 3*n;
edgeNodeMap[pair<int,int>(1,5)] = 7 + 4*n;
edgeNodeMap[pair<int,int>(2,5)] = 7 + 5*n;
edgeNodeMap[pair<int,int>(3,6)] = 7 + 6*n;
edgeNodeMap[pair<int,int>(4,6)] = 7 + 7*n;
edgeNodeMap[pair<int,int>(5,6)] = 7 + 8*n;
// Add vertices
for (int i = 0; i < 6; ++i)
{
m_vertex.push_back(pNodeList[i]);
}
int eid = 0;
// Create edges (with corresponding set of edge points)
for (it = edgeNodeMap.begin(); it != edgeNodeMap.end(); ++it)
{
vector<NodeSharedPtr> edgeNodes;
if (m_conf.m_order > 1) {
for (int j = it->second; j < it->second + n; ++j) {
edgeNodes.push_back(pNodeList[j-1]);
}
}
m_edge.push_back(EdgeSharedPtr(
new Edge(pNodeList[it->first.first-1],
pNodeList[it->first.second-1],
edgeNodes,
m_edge.back()->m_id = eid++;
}
{
}
// Create faces
int face_ids[5][4] = {
{0,1,2,3},{0,1,4,-1},{1,2,5,4},{3,2,5,-1},{0,3,5,4}};
int face_edges[5][4];
int faceoffset = 0;
for (int j = 0; j < 5; ++j)
{
vector<NodeSharedPtr> faceVertices;
vector<EdgeSharedPtr> faceEdges;
vector<NodeSharedPtr> faceNodes;
int nEdge = 3 - (j % 2 - 1);
for (int k = 0; k < nEdge; ++k)
{
faceVertices.push_back(m_vertex[face_ids[j][k]]);
NodeSharedPtr a = m_vertex[face_ids[j][k]];
NodeSharedPtr b = m_vertex[face_ids[j][(k+1) % nEdge]];
for (unsigned int i = 0; i < m_edge.size(); ++i)
{
if ((m_edge[i]->m_n1 == a && m_edge[i]->m_n2 == b) ||
(m_edge[i]->m_n1 == b && m_edge[i]->m_n2 == a))
{
faceEdges.push_back(m_edge[i]);
face_edges[j][k] = i;
break;
}
}
}
{
int facenodes = j%2==0 ? n*n : n*(n-1)/2;
faceoffset += facenodes;
int N = 6 + 9*n + faceoffset;
for (int i = 0; i < facenodes; ++i)
{
faceNodes.push_back(pNodeList[N+i]);
}
}
m_face.push_back(FaceSharedPtr(
new Face(faceVertices, faceNodes, faceEdges, m_conf.m_faceCurveType)));
}
// Re-order edge array to be consistent with Nektar++ ordering.
vector<EdgeSharedPtr> tmp(9);
tmp[0] = m_edge[face_edges[0][0]];
tmp[1] = m_edge[face_edges[0][1]];
tmp[2] = m_edge[face_edges[0][2]];
tmp[3] = m_edge[face_edges[0][3]];
tmp[4] = m_edge[face_edges[1][2]];
tmp[5] = m_edge[face_edges[1][1]];
tmp[6] = m_edge[face_edges[2][1]];
tmp[7] = m_edge[face_edges[3][2]];
tmp[8] = m_edge[face_edges[4][2]];
m_edge = tmp;
}
Nektar::Utilities::Prism::Prism ( const Prism pSrc)
virtual Nektar::Utilities::Prism::~Prism ( )
inlinevirtual

Definition at line 1446 of file MeshElements.h.

{}

Member Function Documentation

void Nektar::Utilities::Prism::Complete ( int  order)
virtual

Reimplemented from Nektar::Utilities::Element.

Definition at line 1496 of file MeshElements.cpp.

References Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eNodalPrismEvenlySpaced, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, GetGeom(), Nektar::Utilities::Element::m_conf, Nektar::Utilities::Element::m_edge, Nektar::Utilities::Element::m_face, Nektar::Utilities::ElmtConfig::m_faceNodes, Nektar::Utilities::ElmtConfig::m_order, Nektar::Utilities::ElmtConfig::m_volumeNodes, and Nektar::Utilities::Element::m_volumeNodes.

{
int i, j, pos;
// Create basis key for a nodal tetrahedron.
LibUtilities::BasisKey B0(
LibUtilities::PointsKey(
LibUtilities::BasisKey B1(
LibUtilities::PointsKey(
LibUtilities::BasisKey B2(
LibUtilities::PointsKey(
// Create a standard nodal prism in order to get the Vandermonde
// matrix to perform interpolation to nodal points.
MemoryManager<StdRegions::StdNodalPrismExp>::AllocateSharedPtr(
Array<OneD, NekDouble> x, y, z;
nodalPrism->GetNodalPoints(x,y,z);
boost::dynamic_pointer_cast<SpatialDomains::PrismGeom>(
this->GetGeom(3));
// Create basis key for a prism.
LibUtilities::BasisKey C0(
LibUtilities::PointsKey(
LibUtilities::BasisKey C1(
LibUtilities::PointsKey(
LibUtilities::BasisKey C2(
LibUtilities::PointsKey(
// Create a prism.
MemoryManager<LocalRegions::PrismExp>::AllocateSharedPtr(
C0, C1, C2, geom);
// Get coordinate array for tetrahedron.
int nqtot = prism->GetTotPoints();
Array<OneD, NekDouble> alloc(6*nqtot);
Array<OneD, NekDouble> xi(alloc);
Array<OneD, NekDouble> yi(alloc+ nqtot);
Array<OneD, NekDouble> zi(alloc+2*nqtot);
Array<OneD, NekDouble> xo(alloc+3*nqtot);
Array<OneD, NekDouble> yo(alloc+4*nqtot);
Array<OneD, NekDouble> zo(alloc+5*nqtot);
Array<OneD, NekDouble> tmp;
prism->GetCoords(xi, yi, zi);
for (i = 0; i < 3; ++i)
{
Array<OneD, NekDouble> coeffs(nodalPrism->GetNcoeffs());
prism->FwdTrans(alloc+i*nqtot, coeffs);
// Apply Vandermonde matrix to project onto nodal space.
nodalPrism->ModalToNodal(coeffs, tmp=alloc+(i+3)*nqtot);
}
// Now extract points from the co-ordinate arrays into the
// edge/face/volume nodes. First, extract edge-interior nodes.
for (i = 0; i < 9; ++i)
{
pos = 6 + i*(order-1);
m_edge[i]->m_edgeNodes.clear();
for (j = 0; j < order-1; ++j)
{
m_edge[i]->m_edgeNodes.push_back(
NodeSharedPtr(new Node(0, xo[pos+j], yo[pos+j], zo[pos+j])));
}
}
// Now extract face-interior nodes.
pos = 6 + 9*(order-1);
for (i = 0; i < 5; ++i)
{
int facesize = i % 2 ? (order-2)*(order-1)/2 : (order-1)*(order-1);
m_face[i]->m_faceNodes.clear();
for (j = 0; j < facesize; ++j)
{
m_face[i]->m_faceNodes.push_back(
NodeSharedPtr(new Node(0, xo[pos+j], yo[pos+j], zo[pos+j])));
}
pos += facesize;
}
// Finally extract volume nodes.
for (i = pos; i < (order+1)*(order+1)*(order+2)/2; ++i)
{
m_volumeNodes.push_back(
NodeSharedPtr(new Node(0, xo[i], yo[i], zo[i])));
}
m_conf.m_order = order;
}
static ElementSharedPtr Nektar::Utilities::Prism::create ( ElmtConfig  pConf,
std::vector< NodeSharedPtr pNodeList,
std::vector< int >  pTagList 
)
inlinestatic

Creates an instance of this class.

Definition at line 1425 of file MeshElements.h.

References Nektar::Utilities::Element::GetFaceList().

{
ElementSharedPtr e = boost::shared_ptr<Element>(
new Prism(pConf, pNodeList, pTagList));
vector<FaceSharedPtr> faces = e->GetFaceList();
for (int i = 0; i < faces.size(); ++i)
{
faces[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e,i));
}
return e;
}
SpatialDomains::GeometrySharedPtr Nektar::Utilities::Prism::GetGeom ( int  coordDim)
virtual

Generate a Nektar++ geometry object for this element.

Reimplemented from Nektar::Utilities::Element.

Definition at line 1477 of file MeshElements.cpp.

References Nektar::Utilities::Element::m_face.

Referenced by Complete().

{
for (int i = 0; i < 5; ++i)
{
faces[i] = m_face[i]->GetGeom(coordDim);
}
ret = MemoryManager<SpatialDomains::PrismGeom>::
AllocateSharedPtr(faces);
return ret;
}
unsigned int Nektar::Utilities::Prism::GetNumNodes ( ElmtConfig  pConf)
static

Return the number of nodes defining a prism.

Definition at line 1466 of file MeshElements.cpp.

References Nektar::Utilities::ElmtConfig::m_faceNodes, Nektar::Utilities::ElmtConfig::m_order, and Nektar::Utilities::ElmtConfig::m_volumeNodes.

Referenced by Nektar::Utilities::InputGmsh::GetNnodes().

{
int n = pConf.m_order;
if (pConf.m_faceNodes && pConf.m_volumeNodes)
return (n+1)*(n+1)*(n+2)/2;
else if (pConf.m_faceNodes && !pConf.m_volumeNodes)
return 3*(n+1)*(n+1)+2*(n+1)*(n+2)/2-9*(n+1)+6;
else
return 9*(n+1)-12;
}
void Nektar::Utilities::Prism::OrientPrism ( )
protected

Orient prism to align degenerate vertices.

Orientation of prismatric elements is required so that the singular vertices of triangular faces (which occur as a part of the collapsed co-ordinate system) align. The algorithm is based on that used in T. Warburton's thesis and in the original Nektar source.

First the points are re-ordered so that the highest global IDs represent the two singular points of the prism. Then, if necessary, the nodes are rotated either clockwise or counter-clockwise (w.r.t to the p-r plane) to correctly align the prism. The #orientation variable is set to:

  • 0 if the prism is not rotated;
  • 1 if the prism is rotated clockwise;
  • 2 if the prism is rotated counter-clockwise.

This is necessary for some input modules (e.g. #InputNek) which add high-order information or bounary conditions to faces.

Definition at line 1627 of file MeshElements.cpp.

References m_orientation, and Nektar::Utilities::Element::m_vertex.

Referenced by Prism().

{
int lid[6], gid[6];
// Re-order vertices.
for (int i = 0; i < 6; ++i)
{
lid[i] = i;
gid[i] = m_vertex[i]->m_id;
}
gid[0] = gid[3] = max(gid[0], gid[3]);
gid[1] = gid[2] = max(gid[1], gid[2]);
gid[4] = gid[5] = max(gid[4], gid[5]);
for (int i = 1; i < 6; ++i)
{
if (gid[0] < gid[i])
{
swap(gid[i], gid[0]);
swap(lid[i], lid[0]);
}
}
if (lid[0] == 4 || lid[0] == 5)
{
}
else if (lid[0] == 1 || lid[0] == 2)
{
// Rotate prism clockwise in p-r plane
vector<NodeSharedPtr> vertexmap(6);
vertexmap[0] = m_vertex[4];
vertexmap[1] = m_vertex[0];
vertexmap[2] = m_vertex[3];
vertexmap[3] = m_vertex[5];
vertexmap[4] = m_vertex[1];
vertexmap[5] = m_vertex[2];
m_vertex = vertexmap;
}
else if (lid[0] == 0 || lid[0] == 3)
{
// Rotate prism counter-clockwise in p-r plane
vector<NodeSharedPtr> vertexmap(6);
vertexmap[0] = m_vertex[1];
vertexmap[1] = m_vertex[4];
vertexmap[2] = m_vertex[5];
vertexmap[3] = m_vertex[2];
vertexmap[4] = m_vertex[0];
vertexmap[5] = m_vertex[3];
m_vertex = vertexmap;
}
else
{
cerr << "Warning: possible prism orientation problem." << endl;
}
}

Member Data Documentation

unsigned int Nektar::Utilities::Prism::m_orientation

Orientation of prism; unchanged = 0; clockwise = 1; counter-clockwise = 2. This is set by OrientPrism.

Definition at line 1457 of file MeshElements.h.

Referenced by OrientPrism().

LibUtilities::ShapeType Nektar::Utilities::Prism::m_type
static
Initial value:
RegisterCreatorFunction(LibUtilities::ePrism, Prism::create, "Prism")

Element type.

Definition at line 1440 of file MeshElements.h.