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::NekMeshUtils::Prism Class Reference

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

#include <Prism.h>

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

Public Member Functions

NEKMESHUTILS_EXPORT Prism (ElmtConfig pConf, std::vector< NodeSharedPtr > pNodeList, std::vector< int > pTagList)
 Create a prism element. More...
 
NEKMESHUTILS_EXPORT Prism (const Prism &pSrc)
 
virtual NEKMESHUTILS_EXPORT ~Prism ()
 
virtual NEKMESHUTILS_EXPORT
SpatialDomains::GeometrySharedPtr 
GetGeom (int coordDim)
 Generate a Nektar++ geometry object for this element. More...
 
virtual NEKMESHUTILS_EXPORT void Complete (int order)
 
- Public Member Functions inherited from Nektar::NekMeshUtils::Element
NEKMESHUTILS_EXPORT Element (ElmtConfig pConf, unsigned int pNumNodes, unsigned int pGotNodes)
 
NEKMESHUTILS_EXPORT unsigned int GetId () const
 Returns the ID of the element (or associated edge or face for boundary elements). More...
 
NEKMESHUTILS_EXPORT unsigned int GetDim () const
 Returns the expansion dimension of the element. More...
 
NEKMESHUTILS_EXPORT ElmtConfig GetConf () const
 Returns the configuration of the element. More...
 
NEKMESHUTILS_EXPORT std::string GetTag () const
 Returns the tag which defines the element shape. More...
 
NEKMESHUTILS_EXPORT NodeSharedPtr GetVertex (unsigned int i) const
 Access a vertex node. More...
 
NEKMESHUTILS_EXPORT EdgeSharedPtr GetEdge (unsigned int i) const
 Access an edge. More...
 
NEKMESHUTILS_EXPORT FaceSharedPtr GetFace (unsigned int i) const
 Access a face. More...
 
NEKMESHUTILS_EXPORT
std::vector< NodeSharedPtr
GetVertexList () const
 Access the list of vertex nodes. More...
 
NEKMESHUTILS_EXPORT
std::vector< EdgeSharedPtr
GetEdgeList () const
 Access the list of edges. More...
 
NEKMESHUTILS_EXPORT
std::vector< FaceSharedPtr
GetFaceList () const
 Access the list of faces. More...
 
NEKMESHUTILS_EXPORT
std::vector< NodeSharedPtr
GetVolumeNodes () const
 Access the list of volume nodes. More...
 
NEKMESHUTILS_EXPORT void SetVolumeNodes (std::vector< NodeSharedPtr > &nodes)
 
NEKMESHUTILS_EXPORT
LibUtilities::PointsType 
GetCurveType () const
 
NEKMESHUTILS_EXPORT void SetCurveType (LibUtilities::PointsType cT)
 
NEKMESHUTILS_EXPORT unsigned int GetNodeCount () const
 Returns the total number of nodes (vertices, edge nodes and face nodes and volume nodes). More...
 
NEKMESHUTILS_EXPORT
std::vector< int > 
GetTagList () const
 Access the list of tags associated with this element. More...
 
NEKMESHUTILS_EXPORT unsigned int GetVertexCount () const
 Returns the number of vertices. More...
 
NEKMESHUTILS_EXPORT unsigned int GetEdgeCount () const
 Returns the number of edges. More...
 
NEKMESHUTILS_EXPORT unsigned int GetFaceCount () const
 Returns the number of faces. More...
 
NEKMESHUTILS_EXPORT void SetId (unsigned int p)
 Change the ID of the element. More...
 
NEKMESHUTILS_EXPORT void SetVertex (unsigned int p, NodeSharedPtr pNew)
 Replace a vertex with another vertex object. More...
 
NEKMESHUTILS_EXPORT void SetEdge (unsigned int p, EdgeSharedPtr pNew)
 Replace an edge with another edge object. More...
 
NEKMESHUTILS_EXPORT void SetFace (unsigned int p, FaceSharedPtr pNew)
 Replace a face with another face object. More...
 
NEKMESHUTILS_EXPORT void SetEdgeLink (EdgeSharedPtr pLink)
 Set a correspondence between this element and an edge (2D boundary element). More...
 
NEKMESHUTILS_EXPORT EdgeSharedPtr GetEdgeLink ()
 Get correspondence between this element and an edge. More...
 
NEKMESHUTILS_EXPORT void SetFaceLink (FaceSharedPtr pLink)
 Set a correspondence between this element and a face (3D boundary element). More...
 
NEKMESHUTILS_EXPORT FaceSharedPtr GetFaceLink ()
 Get correspondence between this element and a face. More...
 
NEKMESHUTILS_EXPORT 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...
 
NEKMESHUTILS_EXPORT int GetBoundaryLink (int i)
 Get the location of the boundary face/edge i for this element. More...
 
NEKMESHUTILS_EXPORT void SetTagList (const std::vector< int > &tags)
 Set the list of tags associated with this element. More...
 
virtual NEKMESHUTILS_EXPORT
std::string 
GetXmlString () const
 Generate a list of vertices (1D), edges (2D), or faces (3D). More...
 
NEKMESHUTILS_EXPORT void GetCurvedNodes (std::vector< NodeSharedPtr > &nodeList) const
 
NEKMESHUTILS_EXPORT std::string GetXmlCurveString () const
 Generates a string listing the coordinates of all nodes associated with this element. More...
 
NEKMESHUTILS_EXPORT int GetMaxOrder ()
 Obtain the order of an element by looking at edges. More...
 
NEKMESHUTILS_EXPORT 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 NEKMESHUTILS_EXPORT
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::NekMeshUtils::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 50 of file Prism.h.

Constructor & Destructor Documentation

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

Create a prism element.

Definition at line 55 of file Prism.cpp.

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

Referenced by create().

58  : Element(pConf, GetNumNodes(pConf), pNodeList.size())
59 {
60  m_tag = "R";
61  m_dim = 3;
62  m_taglist = pTagList;
63  int n = m_conf.m_order - 1;
64 
65  // Create a map to relate edge nodes to a pair of vertices
66  // defining an edge. This is based on the ordering produced by
67  // gmsh.
68  map<pair<int, int>, int> edgeNodeMap;
69  map<pair<int, int>, int>::iterator it;
70 
71  // This edge-node map is based on Nektar++ ordering.
72  edgeNodeMap[pair<int, int>(1, 2)] = 7;
73  edgeNodeMap[pair<int, int>(2, 3)] = 7 + n;
74  edgeNodeMap[pair<int, int>(4, 3)] = 7 + 2 * n;
75  edgeNodeMap[pair<int, int>(1, 4)] = 7 + 3 * n;
76  edgeNodeMap[pair<int, int>(1, 5)] = 7 + 4 * n;
77  edgeNodeMap[pair<int, int>(2, 5)] = 7 + 5 * n;
78  edgeNodeMap[pair<int, int>(3, 6)] = 7 + 6 * n;
79  edgeNodeMap[pair<int, int>(4, 6)] = 7 + 7 * n;
80  edgeNodeMap[pair<int, int>(5, 6)] = 7 + 8 * n;
81 
82  // Add vertices
83  for (int i = 0; i < 6; ++i)
84  {
85  m_vertex.push_back(pNodeList[i]);
86  }
87 
88  int eid = 0;
89  // Create edges (with corresponding set of edge points)
90  for (it = edgeNodeMap.begin(); it != edgeNodeMap.end(); ++it)
91  {
92  vector<NodeSharedPtr> edgeNodes;
93  if (m_conf.m_order > 1)
94  {
95  for (int j = it->second; j < it->second + n; ++j)
96  {
97  edgeNodes.push_back(pNodeList[j - 1]);
98  }
99  }
100  m_edge.push_back(EdgeSharedPtr(new Edge(pNodeList[it->first.first - 1],
101  pNodeList[it->first.second - 1],
102  edgeNodes,
104  m_edge.back()->m_id = eid++;
105  }
106 
107  if (m_conf.m_reorient)
108  {
109  OrientPrism();
110  }
111  else
112  {
113  m_orientation = 0;
114  }
115 
116  // Create faces
117  int face_ids[5][4] = {
118  {0, 1, 2, 3}, {0, 1, 4, -1}, {1, 2, 5, 4}, {3, 2, 5, -1}, {0, 3, 5, 4}};
119  int face_edges[5][4];
120 
121  int face_offset[5];
122  face_offset[0] = 6 + 9 * n;
123  for (int j = 0; j < 4; ++j)
124  {
125  int facenodes = j % 2 == 0 ? n * n : n * (n - 1) / 2;
126  face_offset[j + 1] = face_offset[j] + facenodes;
127  }
128 
129  for (int j = 0; j < 5; ++j)
130  {
131  vector<NodeSharedPtr> faceVertices;
132  vector<EdgeSharedPtr> faceEdges;
133  vector<NodeSharedPtr> faceNodes;
134  int nEdge = 3 - (j % 2 - 1);
135 
136  for (int k = 0; k < nEdge; ++k)
137  {
138  faceVertices.push_back(m_vertex[face_ids[j][k]]);
139  NodeSharedPtr a = m_vertex[face_ids[j][k]];
140  NodeSharedPtr b = m_vertex[face_ids[j][(k + 1) % nEdge]];
141  unsigned int i;
142  for (i = 0; i < m_edge.size(); ++i)
143  {
144  if ((m_edge[i]->m_n1->m_id == a->m_id &&
145  m_edge[i]->m_n2->m_id == b->m_id) ||
146  (m_edge[i]->m_n1->m_id == b->m_id &&
147  m_edge[i]->m_n2->m_id == a->m_id))
148  {
149  faceEdges.push_back(m_edge[i]);
150  face_edges[j][k] = i;
151  break;
152  }
153  }
154 
155  if (i == m_edge.size())
156  {
157  face_edges[j][k] = -1;
158  }
159  }
160 
161  if (m_conf.m_faceNodes)
162  {
163  int face = j, facenodes;
164 
165  if (j % 2 == 0)
166  {
167  facenodes = n * n;
168  if (m_orientation == 1)
169  {
170  face = (face + 4) % 6;
171  }
172  else if (m_orientation == 2)
173  {
174  face = (face + 2) % 6;
175  }
176  }
177  else
178  {
179  // TODO: need to rotate these too.
180  facenodes = n * (n - 1) / 2;
181  }
182 
183  for (int i = 0; i < facenodes; ++i)
184  {
185  faceNodes.push_back(pNodeList[face_offset[face] + i]);
186  }
187  }
188  m_face.push_back(FaceSharedPtr(new Face(
189  faceVertices, faceNodes, faceEdges, m_conf.m_faceCurveType)));
190  }
191 
192  // Re-order edge array to be consistent with Nektar++ ordering.
193  vector<EdgeSharedPtr> tmp(9);
194  ASSERTL1(face_edges[0][0] != -1, "face_edges[0][0] == -1");
195  tmp[0] = m_edge[face_edges[0][0]];
196  ASSERTL1(face_edges[0][1] != -1, "face_edges[0][1] == -1");
197  tmp[1] = m_edge[face_edges[0][1]];
198  ASSERTL1(face_edges[0][2] != -1, "face_edges[0][2] == -1");
199  tmp[2] = m_edge[face_edges[0][2]];
200  ASSERTL1(face_edges[0][3] != -1, "face_edges[0][3] == -1");
201  tmp[3] = m_edge[face_edges[0][3]];
202  ASSERTL1(face_edges[1][2] != -1, "face_edges[1][2] == -1");
203  tmp[4] = m_edge[face_edges[1][2]];
204  ASSERTL1(face_edges[1][1] != -1, "face_edges[1][1] == -1");
205  tmp[5] = m_edge[face_edges[1][1]];
206  ASSERTL1(face_edges[2][1] != -1, "face_edges[2][1] == -1");
207  tmp[6] = m_edge[face_edges[2][1]];
208  ASSERTL1(face_edges[3][2] != -1, "face_edges[3][2] == -1");
209  tmp[7] = m_edge[face_edges[3][2]];
210  ASSERTL1(face_edges[4][2] != -1, "face_edges[4][2] == -1");
211  tmp[8] = m_edge[face_edges[4][2]];
212  m_edge = tmp;
213 }
bool m_faceNodes
Denotes whether the element contains face nodes. For 2D elements, if this is true then the element co...
Definition: Element.h:89
LibUtilities::PointsType m_faceCurveType
Distribution of points in faces.
Definition: Element.h:103
unsigned int m_orientation
Definition: Prism.h:88
ElmtConfig m_conf
Contains configuration of the element.
Definition: Element.h:493
void OrientPrism()
Orient prism to align degenerate vertices.
Definition: Prism.cpp:386
std::vector< int > m_taglist
List of integers specifying properties of the element.
Definition: Element.h:497
LibUtilities::PointsType m_edgeCurveType
Distribution of points in edges.
Definition: Element.h:101
unsigned int m_order
Order of the element.
Definition: Element.h:96
std::vector< NodeSharedPtr > m_vertex
List of element vertex nodes.
Definition: Element.h:499
unsigned int m_dim
Dimension of the element.
Definition: Element.h:491
static NEKMESHUTILS_EXPORT unsigned int GetNumNodes(ElmtConfig pConf)
Return the number of nodes defining a prism.
Definition: Prism.cpp:218
std::vector< EdgeSharedPtr > m_edge
List of element edges.
Definition: Element.h:501
boost::shared_ptr< Node > NodeSharedPtr
Definition: Node.h:50
std::string m_tag
Tag character describing the element.
Definition: Element.h:495
boost::shared_ptr< Edge > EdgeSharedPtr
Shared pointer to an edge.
Definition: Edge.h:196
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
NEKMESHUTILS_EXPORT Element(ElmtConfig pConf, unsigned int pNumNodes, unsigned int pGotNodes)
Definition: Element.cpp:54
std::vector< FaceSharedPtr > m_face
List of element faces.
Definition: Element.h:503
bool m_reorient
Denotes whether the element needs to be re-orientated for a spectral element framework.
Definition: Element.h:99
boost::shared_ptr< Face > FaceSharedPtr
Shared pointer to a face.
Definition: Face.h:378
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
NEKMESHUTILS_EXPORT Nektar::NekMeshUtils::Prism::Prism ( const Prism pSrc)
virtual NEKMESHUTILS_EXPORT Nektar::NekMeshUtils::Prism::~Prism ( )
inlinevirtual

Definition at line 74 of file Prism.h.

75  {
76  }

Member Function Documentation

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

Reimplemented from Nektar::NekMeshUtils::Element.

Definition at line 248 of file Prism.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::NekMeshUtils::Element::m_conf, Nektar::NekMeshUtils::Element::m_edge, Nektar::NekMeshUtils::Element::m_face, Nektar::NekMeshUtils::ElmtConfig::m_faceNodes, Nektar::NekMeshUtils::ElmtConfig::m_order, Nektar::NekMeshUtils::ElmtConfig::m_volumeNodes, and Nektar::NekMeshUtils::Element::m_volumeNodes.

249 {
250  int i, j, pos;
251 
252  // Create basis key for a nodal tetrahedron.
253  LibUtilities::BasisKey B0(
255  order + 1,
256  LibUtilities::PointsKey(order + 1,
258  LibUtilities::BasisKey B1(
260  order + 1,
261  LibUtilities::PointsKey(order + 1,
263  LibUtilities::BasisKey B2(
265  order + 1,
266  LibUtilities::PointsKey(order + 1,
268 
269  // Create a standard nodal prism in order to get the Vandermonde
270  // matrix to perform interpolation to nodal points.
274 
275  Array<OneD, NekDouble> x, y, z;
276  nodalPrism->GetNodalPoints(x, y, z);
277 
279  boost::dynamic_pointer_cast<SpatialDomains::PrismGeom>(
280  this->GetGeom(3));
281 
282  // Create basis key for a prism.
283  LibUtilities::BasisKey C0(
285  order + 1,
286  LibUtilities::PointsKey(order + 1,
288  LibUtilities::BasisKey C1(
290  order + 1,
291  LibUtilities::PointsKey(order + 1,
293  LibUtilities::BasisKey C2(
295  order + 1,
296  LibUtilities::PointsKey(order + 1,
298 
299  // Create a prism.
302  C0, C1, C2, geom);
303 
304  // Get coordinate array for tetrahedron.
305  int nqtot = prism->GetTotPoints();
306  Array<OneD, NekDouble> alloc(6 * nqtot);
307  Array<OneD, NekDouble> xi(alloc);
308  Array<OneD, NekDouble> yi(alloc + nqtot);
309  Array<OneD, NekDouble> zi(alloc + 2 * nqtot);
310  Array<OneD, NekDouble> xo(alloc + 3 * nqtot);
311  Array<OneD, NekDouble> yo(alloc + 4 * nqtot);
312  Array<OneD, NekDouble> zo(alloc + 5 * nqtot);
313  Array<OneD, NekDouble> tmp;
314 
315  prism->GetCoords(xi, yi, zi);
316 
317  for (i = 0; i < 3; ++i)
318  {
319  Array<OneD, NekDouble> coeffs(nodalPrism->GetNcoeffs());
320  prism->FwdTrans(alloc + i * nqtot, coeffs);
321  // Apply Vandermonde matrix to project onto nodal space.
322  nodalPrism->ModalToNodal(coeffs, tmp = alloc + (i + 3) * nqtot);
323  }
324 
325  // Now extract points from the co-ordinate arrays into the
326  // edge/face/volume nodes. First, extract edge-interior nodes.
327  for (i = 0; i < 9; ++i)
328  {
329  pos = 6 + i * (order - 1);
330  m_edge[i]->m_edgeNodes.clear();
331  for (j = 0; j < order - 1; ++j)
332  {
333  m_edge[i]->m_edgeNodes.push_back(NodeSharedPtr(
334  new Node(0, xo[pos + j], yo[pos + j], zo[pos + j])));
335  }
336  }
337 
338  // Now extract face-interior nodes.
339  pos = 6 + 9 * (order - 1);
340  for (i = 0; i < 5; ++i)
341  {
342  int facesize =
343  i % 2 ? (order - 2) * (order - 1) / 2 : (order - 1) * (order - 1);
344  m_face[i]->m_faceNodes.clear();
345  for (j = 0; j < facesize; ++j)
346  {
347  m_face[i]->m_faceNodes.push_back(NodeSharedPtr(
348  new Node(0, xo[pos + j], yo[pos + j], zo[pos + j])));
349  }
350  pos += facesize;
351  }
352 
353  // Finally extract volume nodes.
354  for (i = pos; i < (order + 1) * (order + 1) * (order + 2) / 2; ++i)
355  {
356  m_volumeNodes.push_back(
357  NodeSharedPtr(new Node(0, xo[i], yo[i], zo[i])));
358  }
359 
360  m_conf.m_order = order;
361  m_conf.m_faceNodes = true;
362  m_conf.m_volumeNodes = true;
363 }
bool m_faceNodes
Denotes whether the element contains face nodes. For 2D elements, if this is true then the element co...
Definition: Element.h:89
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
virtual NEKMESHUTILS_EXPORT SpatialDomains::GeometrySharedPtr GetGeom(int coordDim)
Generate a Nektar++ geometry object for this element.
Definition: Prism.cpp:230
ElmtConfig m_conf
Contains configuration of the element.
Definition: Element.h:493
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:57
unsigned int m_order
Order of the element.
Definition: Element.h:96
Principle Orthogonal Functions .
Definition: BasisType.h:47
boost::shared_ptr< StdNodalPrismExp > StdNodalPrismExpSharedPtr
bool m_volumeNodes
Denotes whether the element contains volume (i.e. interior) nodes. These are not supported by either ...
Definition: Element.h:94
std::vector< EdgeSharedPtr > m_edge
List of element edges.
Definition: Element.h:501
boost::shared_ptr< Node > NodeSharedPtr
Definition: Node.h:50
Principle Orthogonal Functions .
Definition: BasisType.h:46
std::vector< NodeSharedPtr > m_volumeNodes
List of element volume nodes.
Definition: Element.h:505
boost::shared_ptr< PrismExp > PrismExpSharedPtr
Definition: PrismExp.h:217
boost::shared_ptr< PrismGeom > PrismGeomSharedPtr
Definition: PrismGeom.h:109
3D Evenly-spaced points on a Prism
Definition: PointsType.h:73
std::vector< FaceSharedPtr > m_face
List of element faces.
Definition: Element.h:503
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50
static ElementSharedPtr Nektar::NekMeshUtils::Prism::create ( ElmtConfig  pConf,
std::vector< NodeSharedPtr pNodeList,
std::vector< int >  pTagList 
)
inlinestatic

Creates an instance of this class.

Definition at line 54 of file Prism.h.

References Prism().

57  {
59  boost::shared_ptr<Element>(new Prism(pConf, pNodeList, pTagList));
60  vector<FaceSharedPtr> faces = e->GetFaceList();
61  for (int i = 0; i < faces.size(); ++i)
62  {
63  faces[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e, i));
64  }
65  return e;
66  }
NEKMESHUTILS_EXPORT Prism(ElmtConfig pConf, std::vector< NodeSharedPtr > pNodeList, std::vector< int > pTagList)
Create a prism element.
Definition: Prism.cpp:55
boost::shared_ptr< Element > ElementSharedPtr
Definition: Edge.h:52
SpatialDomains::GeometrySharedPtr Nektar::NekMeshUtils::Prism::GetGeom ( int  coordDim)
virtual

Generate a Nektar++ geometry object for this element.

Reimplemented from Nektar::NekMeshUtils::Element.

Definition at line 230 of file Prism.cpp.

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

Referenced by Complete().

231 {
234 
235  for (int i = 0; i < 5; ++i)
236  {
237  faces[i] = m_face[i]->GetGeom(coordDim);
238  }
239 
241 
242  return ret;
243 }
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
boost::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry2D.h:59
boost::shared_ptr< PrismGeom > PrismGeomSharedPtr
Definition: PrismGeom.h:109
std::vector< FaceSharedPtr > m_face
List of element faces.
Definition: Element.h:503
unsigned int Nektar::NekMeshUtils::Prism::GetNumNodes ( ElmtConfig  pConf)
static

Return the number of nodes defining a prism.

Definition at line 218 of file Prism.cpp.

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

219 {
220  int n = pConf.m_order;
221  if (pConf.m_faceNodes && pConf.m_volumeNodes)
222  return (n + 1) * (n + 1) * (n + 2) / 2;
223  else if (pConf.m_faceNodes && !pConf.m_volumeNodes)
224  return 3 * (n + 1) * (n + 1) + 2 * (n + 1) * (n + 2) / 2 - 9 * (n + 1) +
225  6;
226  else
227  return 9 * (n + 1) - 12;
228 }
void Nektar::NekMeshUtils::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 386 of file Prism.cpp.

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

Referenced by Prism().

387 {
388  int lid[6], gid[6];
389 
390  // Re-order vertices.
391  for (int i = 0; i < 6; ++i)
392  {
393  lid[i] = i;
394  gid[i] = m_vertex[i]->m_id;
395  }
396 
397  gid[0] = gid[3] = max(gid[0], gid[3]);
398  gid[1] = gid[2] = max(gid[1], gid[2]);
399  gid[4] = gid[5] = max(gid[4], gid[5]);
400 
401  for (int i = 1; i < 6; ++i)
402  {
403  if (gid[0] < gid[i])
404  {
405  swap(gid[i], gid[0]);
406  swap(lid[i], lid[0]);
407  }
408  }
409 
410  if (lid[0] == 4 || lid[0] == 5)
411  {
412  m_orientation = 0;
413  }
414  else if (lid[0] == 1 || lid[0] == 2)
415  {
416  // Rotate prism clockwise in p-r plane
417  vector<NodeSharedPtr> vertexmap(6);
418  vertexmap[0] = m_vertex[4];
419  vertexmap[1] = m_vertex[0];
420  vertexmap[2] = m_vertex[3];
421  vertexmap[3] = m_vertex[5];
422  vertexmap[4] = m_vertex[1];
423  vertexmap[5] = m_vertex[2];
424  m_vertex = vertexmap;
425  m_orientation = 1;
426  }
427  else if (lid[0] == 0 || lid[0] == 3)
428  {
429  // Rotate prism counter-clockwise in p-r plane
430  vector<NodeSharedPtr> vertexmap(6);
431  vertexmap[0] = m_vertex[1];
432  vertexmap[1] = m_vertex[4];
433  vertexmap[2] = m_vertex[5];
434  vertexmap[3] = m_vertex[2];
435  vertexmap[4] = m_vertex[0];
436  vertexmap[5] = m_vertex[3];
437  m_vertex = vertexmap;
438  m_orientation = 2;
439  }
440  else
441  {
442  cerr << "Warning: possible prism orientation problem." << endl;
443  }
444 }
unsigned int m_orientation
Definition: Prism.h:88
std::vector< NodeSharedPtr > m_vertex
List of element vertex nodes.
Definition: Element.h:499

Member Data Documentation

unsigned int Nektar::NekMeshUtils::Prism::m_orientation

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

Definition at line 88 of file Prism.h.

Referenced by OrientPrism(), and Prism().

LibUtilities::ShapeType Nektar::NekMeshUtils::Prism::m_type
static
Initial value:

Element type.

Definition at line 68 of file Prism.h.