Nektar++
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. More...
 
 Prism (const Prism &pSrc)
 
virtual ~Prism ()
 
virtual SpatialDomains::GeometrySharedPtr GetGeom (int coordDim)
 Generate a Nektar++ geometry object for this element. More...
 
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). More...
 
unsigned int GetDim () const
 Returns the expansion dimension of the element. More...
 
ElmtConfig GetConf () const
 Returns the configuration of the element. More...
 
std::string GetTag () const
 Returns the tag which defines the element shape. More...
 
NodeSharedPtr GetVertex (unsigned int i) const
 Access a vertex node. More...
 
EdgeSharedPtr GetEdge (unsigned int i) const
 Access an edge. More...
 
FaceSharedPtr GetFace (unsigned int i) const
 Access a face. More...
 
std::vector< NodeSharedPtrGetVertexList () const
 Access the list of vertex nodes. More...
 
std::vector< EdgeSharedPtrGetEdgeList () const
 Access the list of edges. More...
 
std::vector< FaceSharedPtrGetFaceList () const
 Access the list of faces. More...
 
std::vector< NodeSharedPtrGetVolumeNodes () const
 Access the list of volume nodes. More...
 
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). More...
 
std::vector< int > GetTagList () const
 Access the list of tags associated with this element. More...
 
unsigned int GetVertexCount () const
 Returns the number of vertices. More...
 
unsigned int GetEdgeCount () const
 Returns the number of edges. More...
 
unsigned int GetFaceCount () const
 Returns the number of faces. More...
 
void SetId (unsigned int p)
 Change the ID of the element. More...
 
void SetVertex (unsigned int p, NodeSharedPtr pNew)
 Replace a vertex with another vertex object. More...
 
void SetEdge (unsigned int p, EdgeSharedPtr pNew)
 Replace an edge with another edge object. More...
 
void SetFace (unsigned int p, FaceSharedPtr pNew)
 Replace a face with another face object. More...
 
void SetEdgeLink (EdgeSharedPtr pLink)
 Set a correspondence between this element and an edge (2D boundary element). More...
 
EdgeSharedPtr GetEdgeLink ()
 Get correspondence between this element and an edge. More...
 
void SetFaceLink (FaceSharedPtr pLink)
 Set a correspondence between this element and a face (3D boundary element). More...
 
FaceSharedPtr GetFaceLink ()
 Get correspondence between this element and a face. More...
 
void SetBoundaryLink (int i, int j)
 Set a correspondence between edge or face i and its representative boundary element m->element[expDim-1][j]. More...
 
int GetBoundaryLink (int i)
 Get the location of the boundary face/edge i for this element. More...
 
void SetTagList (const std::vector< int > &tags)
 Set the list of tags associated with this element. More...
 
virtual std::string GetXmlString () const
 Generate a list of vertices (1D), edges (2D), or faces (3D). More...
 
std::string GetXmlCurveString () const
 Generates a string listing the coordinates of all nodes associated with this element. More...
 
int GetMaxOrder ()
 Obtain the order of an element by looking at edges. More...
 
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. More...
 
static unsigned int GetNumNodes (ElmtConfig pConf)
 Return the number of nodes defining a prism. More...
 

Public Attributes

unsigned int m_orientation
 

Static Public Attributes

static LibUtilities::ShapeType m_type
 Element type. More...
 

Protected Member Functions

void OrientPrism ()
 Orient prism to align degenerate vertices. More...
 

Additional Inherited Members

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

Detailed Description

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

Definition at line 1499 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 1453 of file MeshElements.cpp.

References ASSERTL1, 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, m_orientation, Nektar::Utilities::ElmtConfig::m_reorient, Nektar::Utilities::Element::m_tag, Nektar::Utilities::Element::m_taglist, Nektar::Utilities::Element::m_vertex, and OrientPrism().

1456  : Element(pConf, GetNumNodes(pConf), pNodeList.size())
1457  {
1458  m_tag = "R";
1459  m_dim = 3;
1460  m_taglist = pTagList;
1461  int n = m_conf.m_order-1;
1462 
1463  // Create a map to relate edge nodes to a pair of vertices
1464  // defining an edge. This is based on the ordering produced by
1465  // gmsh.
1466  map<pair<int,int>, int> edgeNodeMap;
1467  map<pair<int,int>, int>::iterator it;
1468 
1469  // This edge-node map is based on Nektar++ ordering.
1470  edgeNodeMap[pair<int,int>(1,2)] = 7;
1471  edgeNodeMap[pair<int,int>(2,3)] = 7 + n;
1472  edgeNodeMap[pair<int,int>(4,3)] = 7 + 2*n;
1473  edgeNodeMap[pair<int,int>(1,4)] = 7 + 3*n;
1474  edgeNodeMap[pair<int,int>(1,5)] = 7 + 4*n;
1475  edgeNodeMap[pair<int,int>(2,5)] = 7 + 5*n;
1476  edgeNodeMap[pair<int,int>(3,6)] = 7 + 6*n;
1477  edgeNodeMap[pair<int,int>(4,6)] = 7 + 7*n;
1478  edgeNodeMap[pair<int,int>(5,6)] = 7 + 8*n;
1479 
1480  // Add vertices
1481  for (int i = 0; i < 6; ++i)
1482  {
1483  m_vertex.push_back(pNodeList[i]);
1484  }
1485 
1486  int eid = 0;
1487  // Create edges (with corresponding set of edge points)
1488  for (it = edgeNodeMap.begin(); it != edgeNodeMap.end(); ++it)
1489  {
1490  vector<NodeSharedPtr> edgeNodes;
1491  if (m_conf.m_order > 1) {
1492  for (int j = it->second; j < it->second + n; ++j) {
1493  edgeNodes.push_back(pNodeList[j-1]);
1494  }
1495  }
1496  m_edge.push_back(EdgeSharedPtr(
1497  new Edge(pNodeList[it->first.first-1],
1498  pNodeList[it->first.second-1],
1499  edgeNodes,
1501  m_edge.back()->m_id = eid++;
1502  }
1503 
1504  if (m_conf.m_reorient)
1505  {
1506  OrientPrism();
1507  }
1508  else
1509  {
1510  m_orientation = 0;
1511  }
1512 
1513  // Create faces
1514  int face_ids[5][4] = {
1515  {0,1,2,3},{0,1,4,-1},{1,2,5,4},{3,2,5,-1},{0,3,5,4}};
1516  int face_edges[5][4];
1517 
1518  int face_offset[5];
1519  face_offset[0] = 6 + 9*n;
1520  for (int j = 0; j < 4; ++j)
1521  {
1522  int facenodes = j % 2 == 0 ? n*n : n*(n-1)/2;
1523  face_offset[j+1] = face_offset[j] + facenodes;
1524  }
1525 
1526  for (int j = 0; j < 5; ++j)
1527  {
1528  vector<NodeSharedPtr> faceVertices;
1529  vector<EdgeSharedPtr> faceEdges;
1530  vector<NodeSharedPtr> faceNodes;
1531  int nEdge = 3 - (j % 2 - 1);
1532 
1533  for (int k = 0; k < nEdge; ++k)
1534  {
1535  faceVertices.push_back(m_vertex[face_ids[j][k]]);
1536  NodeSharedPtr a = m_vertex[face_ids[j][k]];
1537  NodeSharedPtr b = m_vertex[face_ids[j][(k+1) % nEdge]];
1538  unsigned int i;
1539  for (i = 0; i < m_edge.size(); ++i)
1540  {
1541  if ((m_edge[i]->m_n1->m_id == a->m_id
1542  && m_edge[i]->m_n2->m_id == b->m_id) ||
1543  (m_edge[i]->m_n1->m_id == b->m_id
1544  && m_edge[i]->m_n2->m_id == a->m_id))
1545  {
1546  faceEdges.push_back(m_edge[i]);
1547  face_edges[j][k] = i;
1548  break;
1549  }
1550  }
1551 
1552  if(i == m_edge.size())
1553  {
1554  face_edges[j][k] = -1;
1555  }
1556  }
1557 
1558  if (m_conf.m_faceNodes)
1559  {
1560  int face = j, facenodes;
1561 
1562  if (j % 2 == 0)
1563  {
1564  facenodes = n*n;
1565  if (m_orientation == 1)
1566  {
1567  face = (face+4) % 6;
1568  }
1569  else if (m_orientation == 2)
1570  {
1571  face = (face+2) % 6;
1572  }
1573  }
1574  else
1575  {
1576  // TODO: need to rotate these too.
1577  facenodes = n*(n-1)/2;
1578  }
1579 
1580  for (int i = 0; i < facenodes; ++i)
1581  {
1582  faceNodes.push_back(pNodeList[face_offset[face]+i]);
1583  }
1584  }
1585  m_face.push_back(FaceSharedPtr(
1586  new Face(faceVertices, faceNodes, faceEdges, m_conf.m_faceCurveType)));
1587  }
1588 
1589  // Re-order edge array to be consistent with Nektar++ ordering.
1590  vector<EdgeSharedPtr> tmp(9);
1591  ASSERTL1(face_edges[0][0] != -1,"face_edges[0][0] == -1");
1592  tmp[0] = m_edge[face_edges[0][0]];
1593  ASSERTL1(face_edges[0][1] != -1,"face_edges[0][1] == -1");
1594  tmp[1] = m_edge[face_edges[0][1]];
1595  ASSERTL1(face_edges[0][2] != -1,"face_edges[0][2] == -1");
1596  tmp[2] = m_edge[face_edges[0][2]];
1597  ASSERTL1(face_edges[0][3] != -1,"face_edges[0][3] == -1");
1598  tmp[3] = m_edge[face_edges[0][3]];
1599  ASSERTL1(face_edges[1][2] != -1,"face_edges[1][2] == -1");
1600  tmp[4] = m_edge[face_edges[1][2]];
1601  ASSERTL1(face_edges[1][1] != -1,"face_edges[1][1] == -1");
1602  tmp[5] = m_edge[face_edges[1][1]];
1603  ASSERTL1(face_edges[2][1] != -1,"face_edges[2][1] == -1");
1604  tmp[6] = m_edge[face_edges[2][1]];
1605  ASSERTL1(face_edges[3][2] != -1,"face_edges[3][2] == -1");
1606  tmp[7] = m_edge[face_edges[3][2]];
1607  ASSERTL1(face_edges[4][2] != -1,"face_edges[4][2] == -1");
1608  tmp[8] = m_edge[face_edges[4][2]];
1609  m_edge = tmp;
1610  }
static unsigned int GetNumNodes(ElmtConfig pConf)
Return the number of nodes defining a prism.
std::vector< int > m_taglist
List of integers specifying properties of the element.
Definition: MeshElements.h:984
boost::shared_ptr< Edge > EdgeSharedPtr
Shared pointer to an edge.
Definition: MeshElements.h:318
std::vector< FaceSharedPtr > m_face
List of element faces.
Definition: MeshElements.h:990
boost::shared_ptr< Face > FaceSharedPtr
Shared pointer to a face.
Definition: MeshElements.h:550
boost::shared_ptr< Node > NodeSharedPtr
Shared pointer to a Node.
Definition: MeshElements.h:195
bool m_reorient
Denotes whether the element needs to be re-orientated for a spectral element framework.
Definition: MeshElements.h:632
unsigned int m_order
Order of the element.
Definition: MeshElements.h:629
LibUtilities::PointsType m_faceCurveType
Distribution of points in faces.
Definition: MeshElements.h:636
bool m_faceNodes
Denotes whether the element contains face nodes. For 2D elements, if this is true then the element co...
Definition: MeshElements.h:622
std::string m_tag
Tag character describing the element.
Definition: MeshElements.h:982
unsigned int m_dim
Dimension of the element.
Definition: MeshElements.h:978
LibUtilities::PointsType m_edgeCurveType
Distribution of points in edges.
Definition: MeshElements.h:634
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
void OrientPrism()
Orient prism to align degenerate vertices.
std::vector< NodeSharedPtr > m_vertex
List of element vertex nodes.
Definition: MeshElements.h:986
std::vector< EdgeSharedPtr > m_edge
List of element edges.
Definition: MeshElements.h:988
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
Element(ElmtConfig pConf, unsigned int pNumNodes, unsigned int pGotNodes)
ElmtConfig m_conf
Contains configuration of the element.
Definition: MeshElements.h:980
Nektar::Utilities::Prism::Prism ( const Prism pSrc)
virtual Nektar::Utilities::Prism::~Prism ( )
inlinevirtual

Definition at line 1523 of file MeshElements.h.

1523 {}

Member Function Documentation

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

Reimplemented from Nektar::Utilities::Element.

Definition at line 1645 of file MeshElements.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), 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.

1646  {
1647  int i, j, pos;
1648 
1649  // Create basis key for a nodal tetrahedron.
1650  LibUtilities::BasisKey B0(
1651  LibUtilities::eOrtho_A, order+1,
1652  LibUtilities::PointsKey(
1654  LibUtilities::BasisKey B1(
1655  LibUtilities::eOrtho_A, order+1,
1656  LibUtilities::PointsKey(
1658  LibUtilities::BasisKey B2(
1659  LibUtilities::eOrtho_B, order+1,
1660  LibUtilities::PointsKey(
1662 
1663  // Create a standard nodal prism in order to get the Vandermonde
1664  // matrix to perform interpolation to nodal points.
1668 
1669  Array<OneD, NekDouble> x, y, z;
1670  nodalPrism->GetNodalPoints(x,y,z);
1671 
1673  boost::dynamic_pointer_cast<SpatialDomains::PrismGeom>(
1674  this->GetGeom(3));
1675 
1676  // Create basis key for a prism.
1677  LibUtilities::BasisKey C0(
1678  LibUtilities::eOrtho_A, order+1,
1679  LibUtilities::PointsKey(
1681  LibUtilities::BasisKey C1(
1682  LibUtilities::eOrtho_A, order+1,
1683  LibUtilities::PointsKey(
1685  LibUtilities::BasisKey C2(
1686  LibUtilities::eOrtho_B, order+1,
1687  LibUtilities::PointsKey(
1689 
1690  // Create a prism.
1693  C0, C1, C2, geom);
1694 
1695  // Get coordinate array for tetrahedron.
1696  int nqtot = prism->GetTotPoints();
1697  Array<OneD, NekDouble> alloc(6*nqtot);
1698  Array<OneD, NekDouble> xi(alloc);
1699  Array<OneD, NekDouble> yi(alloc+ nqtot);
1700  Array<OneD, NekDouble> zi(alloc+2*nqtot);
1701  Array<OneD, NekDouble> xo(alloc+3*nqtot);
1702  Array<OneD, NekDouble> yo(alloc+4*nqtot);
1703  Array<OneD, NekDouble> zo(alloc+5*nqtot);
1704  Array<OneD, NekDouble> tmp;
1705 
1706  prism->GetCoords(xi, yi, zi);
1707 
1708  for (i = 0; i < 3; ++i)
1709  {
1710  Array<OneD, NekDouble> coeffs(nodalPrism->GetNcoeffs());
1711  prism->FwdTrans(alloc+i*nqtot, coeffs);
1712  // Apply Vandermonde matrix to project onto nodal space.
1713  nodalPrism->ModalToNodal(coeffs, tmp=alloc+(i+3)*nqtot);
1714  }
1715 
1716  // Now extract points from the co-ordinate arrays into the
1717  // edge/face/volume nodes. First, extract edge-interior nodes.
1718  for (i = 0; i < 9; ++i)
1719  {
1720  pos = 6 + i*(order-1);
1721  m_edge[i]->m_edgeNodes.clear();
1722  for (j = 0; j < order-1; ++j)
1723  {
1724  m_edge[i]->m_edgeNodes.push_back(
1725  NodeSharedPtr(new Node(0, xo[pos+j], yo[pos+j], zo[pos+j])));
1726  }
1727  }
1728 
1729  // Now extract face-interior nodes.
1730  pos = 6 + 9*(order-1);
1731  for (i = 0; i < 5; ++i)
1732  {
1733  int facesize = i % 2 ? (order-2)*(order-1)/2 : (order-1)*(order-1);
1734  m_face[i]->m_faceNodes.clear();
1735  for (j = 0; j < facesize; ++j)
1736  {
1737  m_face[i]->m_faceNodes.push_back(
1738  NodeSharedPtr(new Node(0, xo[pos+j], yo[pos+j], zo[pos+j])));
1739  }
1740  pos += facesize;
1741  }
1742 
1743  // Finally extract volume nodes.
1744  for (i = pos; i < (order+1)*(order+1)*(order+2)/2; ++i)
1745  {
1746  m_volumeNodes.push_back(
1747  NodeSharedPtr(new Node(0, xo[i], yo[i], zo[i])));
1748  }
1749 
1750  m_conf.m_order = order;
1751  m_conf.m_faceNodes = true;
1752  m_conf.m_volumeNodes = true;
1753  }
bool m_volumeNodes
Denotes whether the element contains volume (i.e. interior) nodes. These are not supported by either ...
Definition: MeshElements.h:627
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
std::vector< NodeSharedPtr > m_volumeNodes
List of element volume nodes.
Definition: MeshElements.h:992
std::vector< FaceSharedPtr > m_face
List of element faces.
Definition: MeshElements.h:990
boost::shared_ptr< Node > NodeSharedPtr
Shared pointer to a Node.
Definition: MeshElements.h:195
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:57
Principle Orthogonal Functions .
Definition: BasisType.h:47
boost::shared_ptr< StdNodalPrismExp > StdNodalPrismExpSharedPtr
unsigned int m_order
Order of the element.
Definition: MeshElements.h:629
bool m_faceNodes
Denotes whether the element contains face nodes. For 2D elements, if this is true then the element co...
Definition: MeshElements.h:622
Principle Orthogonal Functions .
Definition: BasisType.h:46
boost::shared_ptr< PrismExp > PrismExpSharedPtr
Definition: PrismExp.h:226
virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim)
Generate a Nektar++ geometry object for this element.
boost::shared_ptr< PrismGeom > PrismGeomSharedPtr
Definition: PrismGeom.h:109
3D Evenly-spaced points on a Prism
Definition: PointsType.h:73
std::vector< EdgeSharedPtr > m_edge
List of element edges.
Definition: MeshElements.h:988
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50
ElmtConfig m_conf
Contains configuration of the element.
Definition: MeshElements.h:980
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 1502 of file MeshElements.h.

1506  {
1507  ElementSharedPtr e = boost::shared_ptr<Element>(
1508  new Prism(pConf, pNodeList, pTagList));
1509  vector<FaceSharedPtr> faces = e->GetFaceList();
1510  for (int i = 0; i < faces.size(); ++i)
1511  {
1512  faces[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e,i));
1513  }
1514  return e;
1515  }
boost::shared_ptr< Element > ElementSharedPtr
Shared pointer to an element.
Definition: MeshElements.h:63
Prism(ElmtConfig pConf, std::vector< NodeSharedPtr > pNodeList, std::vector< int > pTagList)
Create a prism element.
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 1626 of file MeshElements.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and Nektar::Utilities::Element::m_face.

Referenced by Complete().

1627  {
1630 
1631  for (int i = 0; i < 5; ++i)
1632  {
1633  faces[i] = m_face[i]->GetGeom(coordDim);
1634  }
1635 
1637  AllocateSharedPtr(faces);
1638 
1639  return ret;
1640  }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
std::vector< FaceSharedPtr > m_face
List of element faces.
Definition: MeshElements.h:990
boost::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry2D.h:59
boost::shared_ptr< PrismGeom > PrismGeomSharedPtr
Definition: PrismGeom.h:109
unsigned int Nektar::Utilities::Prism::GetNumNodes ( ElmtConfig  pConf)
static

Return the number of nodes defining a prism.

Definition at line 1615 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().

1616  {
1617  int n = pConf.m_order;
1618  if (pConf.m_faceNodes && pConf.m_volumeNodes)
1619  return (n+1)*(n+1)*(n+2)/2;
1620  else if (pConf.m_faceNodes && !pConf.m_volumeNodes)
1621  return 3*(n+1)*(n+1)+2*(n+1)*(n+2)/2-9*(n+1)+6;
1622  else
1623  return 9*(n+1)-12;
1624  }
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 1776 of file MeshElements.cpp.

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

Referenced by Prism().

1777  {
1778  int lid[6], gid[6];
1779 
1780  // Re-order vertices.
1781  for (int i = 0; i < 6; ++i)
1782  {
1783  lid[i] = i;
1784  gid[i] = m_vertex[i]->m_id;
1785  }
1786 
1787  gid[0] = gid[3] = max(gid[0], gid[3]);
1788  gid[1] = gid[2] = max(gid[1], gid[2]);
1789  gid[4] = gid[5] = max(gid[4], gid[5]);
1790 
1791  for (int i = 1; i < 6; ++i)
1792  {
1793  if (gid[0] < gid[i])
1794  {
1795  swap(gid[i], gid[0]);
1796  swap(lid[i], lid[0]);
1797  }
1798  }
1799 
1800  if (lid[0] == 4 || lid[0] == 5)
1801  {
1802  m_orientation = 0;
1803  }
1804  else if (lid[0] == 1 || lid[0] == 2)
1805  {
1806  // Rotate prism clockwise in p-r plane
1807  vector<NodeSharedPtr> vertexmap(6);
1808  vertexmap[0] = m_vertex[4];
1809  vertexmap[1] = m_vertex[0];
1810  vertexmap[2] = m_vertex[3];
1811  vertexmap[3] = m_vertex[5];
1812  vertexmap[4] = m_vertex[1];
1813  vertexmap[5] = m_vertex[2];
1814  m_vertex = vertexmap;
1815  m_orientation = 1;
1816  }
1817  else if (lid[0] == 0 || lid[0] == 3)
1818  {
1819  // Rotate prism counter-clockwise in p-r plane
1820  vector<NodeSharedPtr> vertexmap(6);
1821  vertexmap[0] = m_vertex[1];
1822  vertexmap[1] = m_vertex[4];
1823  vertexmap[2] = m_vertex[5];
1824  vertexmap[3] = m_vertex[2];
1825  vertexmap[4] = m_vertex[0];
1826  vertexmap[5] = m_vertex[3];
1827  m_vertex = vertexmap;
1828  m_orientation = 2;
1829  }
1830  else
1831  {
1832  cerr << "Warning: possible prism orientation problem." << endl;
1833  }
1834  }
std::vector< NodeSharedPtr > m_vertex
List of element vertex nodes.
Definition: MeshElements.h:986

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 1534 of file MeshElements.h.

Referenced by OrientPrism(), and Prism().

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

Element type.

Definition at line 1517 of file MeshElements.h.