Nektar++
Public Member Functions | Static Public Member Functions | Public Attributes | Static Public Attributes | Protected Member Functions | Static Private Attributes | 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:
[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 GetCurvedNodes (std::vector< NodeSharedPtr > &nodeList) const
 get list of volume interior nodes More...
 
virtual NEKMESHUTILS_EXPORT StdRegions::Orientation GetEdgeOrient (int edgeId, EdgeSharedPtr edge)
 Get the edge orientation of edge with respect to the local element, which lies at edge index edgeId. More...
 
virtual NEKMESHUTILS_EXPORT void MakeOrder (int order, SpatialDomains::GeometrySharedPtr geom, LibUtilities::PointsType pType, int coordDim, int &id, bool justConfig=false)
 Insert interior (i.e. volume) points into this element to make the geometry an order order representation. More...
 
virtual NEKMESHUTILS_EXPORT int GetFaceVertex (int i, int j)
 Returns the local index of vertex j of face i. More...
 
- 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 LibUtilities::ShapeType GetShapeType () const
 returns the shapetype 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< NodeSharedPtrGetVertexList () const
 Access the list of vertex nodes. More...
 
NEKMESHUTILS_EXPORT std::vector< EdgeSharedPtrGetEdgeList () const
 Access the list of edges. More...
 
NEKMESHUTILS_EXPORT std::vector< FaceSharedPtrGetFaceList () const
 Access the list of faces. More...
 
NEKMESHUTILS_EXPORT std::vector< NodeSharedPtrGetVolumeNodes () 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 ()
 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, bool descend=true)
 Replace a vertex in the element. More...
 
NEKMESHUTILS_EXPORT void SetEdge (unsigned int p, EdgeSharedPtr pNew, bool descend=true)
 Replace an edge in the element. More...
 
NEKMESHUTILS_EXPORT void SetFace (unsigned int p, FaceSharedPtr pNew)
 Replace a face in the element. 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 bool HasBoundaryLinks ()
 Is this element connected to a boundary. 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 ()
 Generate a list of vertices (1D), edges (2D), or faces (3D). More...
 
NEKMESHUTILS_EXPORT std::string GetXmlCurveString ()
 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 bool IsDeformed ()
 Determines whether an element is deformed by inspecting whether there are any edge, face or volume interior nodes. More...
 
NEKMESHUTILS_EXPORT std::pair< Node, NodeGetBoundingBox ()
 Returns the approximate bounding box of this element based on the coordinates of all vertices, edges and faces of the element. Note that this does not robustly take into account the curvature of the element. More...
 
NEKMESHUTILS_EXPORT void Print ()
 
virtual NEKMESHUTILS_EXPORT Array< OneD, NekDoubleNormal (bool inward=false)
 returns the normal to the element More...
 

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
 
- Public Attributes inherited from Nektar::NekMeshUtils::Element
CADObjectSharedPtr m_parentCAD
 

Static Public Attributes

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

Protected Member Functions

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

Static Private Attributes

static int m_faceIds [5][4]
 Vertex IDs that make up prism faces. 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 49 of file Prism.h.

Constructor & Destructor Documentation

◆ Prism() [1/2]

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

Create a prism element.

Definition at line 63 of file Prism.cpp.

References ASSERTL1, Nektar::LibUtilities::eNodalTriEvenlySpaced, Nektar::LibUtilities::ePolyEvenlySpaced, 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, m_faceIds, 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().

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

◆ Prism() [2/2]

NEKMESHUTILS_EXPORT Nektar::NekMeshUtils::Prism::Prism ( const Prism pSrc)

◆ ~Prism()

virtual NEKMESHUTILS_EXPORT Nektar::NekMeshUtils::Prism::~Prism ( )
inlinevirtual

Definition at line 66 of file Prism.h.

References GetCurvedNodes(), GetEdgeOrient(), GetGeom(), GetNumNodes(), MakeOrder(), and NEKMESHUTILS_EXPORT.

67  {
68  }

Member Function Documentation

◆ create()

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 53 of file Prism.h.

56  {
57  return std::make_shared<Prism>(pConf, pNodeList, pTagList);
58  }

◆ GetCurvedNodes()

void Nektar::NekMeshUtils::Prism::GetCurvedNodes ( std::vector< NodeSharedPtr > &  nodeList) const
virtual

get list of volume interior nodes

Reimplemented from Nektar::NekMeshUtils::Element.

Definition at line 354 of file Prism.cpp.

References Nektar::NekMeshUtils::HOTriangle< T >::Align(), Nektar::NekMeshUtils::HOQuadrilateral< T >::Align(), CellMLToNektar.pycml::copy(), Nektar::NekMeshUtils::Element::m_edge, Nektar::NekMeshUtils::Element::m_face, Nektar::NekMeshUtils::Element::m_id, Nektar::NekMeshUtils::Element::m_vertex, Nektar::NekMeshUtils::Element::m_volumeNodes, Nektar::NekMeshUtils::HOTriangle< T >::surfVerts, and Nektar::NekMeshUtils::HOQuadrilateral< T >::surfVerts.

Referenced by ~Prism().

355 {
356  int n = m_edge[0]->GetNodeCount();
357  nodeList.resize(n*n*(n+1)/2);
358 
359  nodeList[0] = m_vertex[0];
360  nodeList[1] = m_vertex[1];
361  nodeList[2] = m_vertex[2];
362  nodeList[3] = m_vertex[3];
363  nodeList[4] = m_vertex[4];
364  nodeList[5] = m_vertex[5];
365  int k = 6;
366 
367  for(int i = 0; i < 4; i++)
368  {
369  bool reverseEdge = m_edge[i]->m_n1 == m_vertex[i];
370  if (reverseEdge)
371  {
372  for(int j = 0; j < n-2; j++)
373  {
374  nodeList[k++] = m_edge[i]->m_edgeNodes[j];
375  }
376  }
377  else
378  {
379  for(int j = n-3; j >= 0; j--)
380  {
381  nodeList[k++] = m_edge[i]->m_edgeNodes[j];
382  }
383  }
384  }
385 
386  for(int i = 4; i < 8; i++)
387  {
388  bool reverseEdge = m_edge[i]->m_n1 == m_vertex[i-4];
389  if (reverseEdge)
390  {
391  for(int j = 0; j < n-2; j++)
392  {
393  nodeList[k++] = m_edge[i]->m_edgeNodes[j];
394  }
395  }
396  else
397  {
398  for(int j = n-3; j >= 0; j--)
399  {
400  nodeList[k++] = m_edge[i]->m_edgeNodes[j];
401  }
402  }
403  }
404  bool reverseEdge = m_edge[8]->m_n1 == m_vertex[4];
405  if (reverseEdge)
406  {
407  for(int j = 0; j < n-2; j++)
408  {
409  nodeList[k++] = m_edge[8]->m_edgeNodes[j];
410  }
411  }
412  else
413  {
414  for(int j = n-3; j >= 0; j--)
415  {
416  nodeList[k++] = m_edge[8]->m_edgeNodes[j];
417  }
418  }
419 
420  vector<vector<int> > ts;
421  {
422  vector<int> t(4);
423  t[0] = m_vertex[0]->m_id;
424  t[1] = m_vertex[1]->m_id;
425  t[2] = m_vertex[2]->m_id;
426  t[3] = m_vertex[3]->m_id;
427  ts.push_back(t);
428  }
429  {
430  vector<int> t(3);
431  t[0] = m_vertex[0]->m_id;
432  t[1] = m_vertex[1]->m_id;
433  t[2] = m_vertex[4]->m_id;
434  ts.push_back(t);
435  }
436  {
437  vector<int> t(4);
438  t[0] = m_vertex[1]->m_id;
439  t[1] = m_vertex[2]->m_id;
440  t[2] = m_vertex[5]->m_id;
441  t[3] = m_vertex[4]->m_id;
442  ts.push_back(t);
443  }
444  {
445  vector<int> t(3);
446  t[0] = m_vertex[3]->m_id;
447  t[1] = m_vertex[2]->m_id;
448  t[2] = m_vertex[5]->m_id;
449  ts.push_back(t);
450  }
451  {
452  vector<int> t(4);
453  t[0] = m_vertex[0]->m_id;
454  t[1] = m_vertex[3]->m_id;
455  t[2] = m_vertex[5]->m_id;
456  t[3] = m_vertex[4]->m_id;
457  ts.push_back(t);
458  }
459 
460  for(int i = 0; i < ts.size(); i++)
461  {
462  if(ts[i].size() == 3)
463  {
464  vector<int> fcid;
465  fcid.push_back(m_face[i]->m_vertexList[0]->m_id);
466  fcid.push_back(m_face[i]->m_vertexList[1]->m_id);
467  fcid.push_back(m_face[i]->m_vertexList[2]->m_id);
468 
469  HOTriangle<NodeSharedPtr> hot(fcid, m_face[i]->m_faceNodes);
470 
471  hot.Align(ts[i]);
472 
473  std::copy(hot.surfVerts.begin(),
474  hot.surfVerts.end(),
475  nodeList.begin() + k);
476  k+= hot.surfVerts.size();
477  }
478  else
479  {
480  vector<int> fcid;
481  fcid.push_back(m_face[i]->m_vertexList[0]->m_id);
482  fcid.push_back(m_face[i]->m_vertexList[1]->m_id);
483  fcid.push_back(m_face[i]->m_vertexList[2]->m_id);
484  fcid.push_back(m_face[i]->m_vertexList[3]->m_id);
485 
486  HOQuadrilateral<NodeSharedPtr> hoq(fcid, m_face[i]->m_faceNodes);
487 
488  hoq.Align(ts[i]);
489 
490  std::copy(hoq.surfVerts.begin(),
491  hoq.surfVerts.end(),
492  nodeList.begin() + k);
493  k+= hoq.surfVerts.size();
494  }
495  }
496 
497  std::copy(m_volumeNodes.begin(),
498  m_volumeNodes.end(),
499  nodeList.begin() + k);
500 }
def copy(self)
Definition: pycml.py:2663
std::vector< NodeSharedPtr > m_vertex
List of element vertex nodes.
Definition: Element.h:468
std::vector< EdgeSharedPtr > m_edge
List of element edges.
Definition: Element.h:470
std::vector< NodeSharedPtr > m_volumeNodes
List of element volume nodes.
Definition: Element.h:474
unsigned int m_id
ID of the element.
Definition: Element.h:458
std::vector< FaceSharedPtr > m_face
List of element faces.
Definition: Element.h:472

◆ GetEdgeOrient()

StdRegions::Orientation Nektar::NekMeshUtils::Prism::GetEdgeOrient ( int  edgeId,
EdgeSharedPtr  edge 
)
virtual

Get the edge orientation of edge with respect to the local element, which lies at edge index edgeId.

Reimplemented from Nektar::NekMeshUtils::Element.

Definition at line 262 of file Prism.cpp.

References ASSERTL1, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eForwards, Nektar::StdRegions::eNoOrientation, and Nektar::NekMeshUtils::Element::m_vertex.

Referenced by ~Prism().

264 {
265  static int edgeVerts[9][2] = {
266  {0,1}, {1,2}, {3,2}, {0,3}, {0,4}, {1,4}, {2,5}, {3,5}, {4,5}
267  };
268 
269  if (edge->m_n1 == m_vertex[edgeVerts[edgeId][0]])
270  {
271  return StdRegions::eForwards;
272  }
273  else if (edge->m_n1 == m_vertex[edgeVerts[edgeId][1]])
274  {
275  return StdRegions::eBackwards;
276  }
277  else
278  {
279  ASSERTL1(false, "Edge is not connected to this quadrilateral.");
280  }
281 
283 }
std::vector< NodeSharedPtr > m_vertex
List of element vertex nodes.
Definition: Element.h:468
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250

◆ GetFaceVertex()

virtual NEKMESHUTILS_EXPORT int Nektar::NekMeshUtils::Prism::GetFaceVertex ( int  i,
int  j 
)
inlinevirtual

Returns the local index of vertex j of face i.

Reimplemented from Nektar::NekMeshUtils::Element.

Definition at line 85 of file Prism.h.

References m_faceIds.

86  {
87  return m_faceIds[i][j];
88  }
static int m_faceIds[5][4]
Vertex IDs that make up prism faces.
Definition: Prism.h:100

◆ GetGeom()

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 245 of file Prism.cpp.

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

Referenced by ~Prism().

246 {
249 
250  for (int i = 0; i < 5; ++i)
251  {
252  faces[i] = m_face[i]->GetGeom(coordDim);
253  }
254 
256  m_id, faces);
257 
258  ret->Setup();
259  return ret;
260 }
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry.h:65
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
unsigned int m_id
ID of the element.
Definition: Element.h:458
std::vector< FaceSharedPtr > m_face
List of element faces.
Definition: Element.h:472
std::shared_ptr< PrismGeom > PrismGeomSharedPtr
Definition: PrismGeom.h:88

◆ GetNumNodes()

unsigned int Nektar::NekMeshUtils::Prism::GetNumNodes ( ElmtConfig  pConf)
static

Return the number of nodes defining a prism.

Definition at line 233 of file Prism.cpp.

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

Referenced by Nektar::Utilities::InputGmsh::GetNnodes(), and ~Prism().

234 {
235  int n = pConf.m_order;
236  if (pConf.m_faceNodes && pConf.m_volumeNodes)
237  return (n + 1) * (n + 1) * (n + 2) / 2;
238  else if (pConf.m_faceNodes && !pConf.m_volumeNodes)
239  return 3 * (n + 1) * (n + 1) + 2 * (n + 1) * (n + 2) / 2 - 9 * (n + 1) +
240  6;
241  else
242  return 9 * (n + 1) - 12;
243 }

◆ MakeOrder()

void Nektar::NekMeshUtils::Prism::MakeOrder ( int  order,
SpatialDomains::GeometrySharedPtr  geom,
LibUtilities::PointsType  edgeType,
int  coordDim,
int &  id,
bool  justConfig = false 
)
virtual

Insert interior (i.e. volume) points into this element to make the geometry an order order representation.

Parameters
orderThe desired polynomial order.
geomThe geometry object used to describe the curvature mapping.
edgeTypeThe points distribution to use on the volume.
coordDimThe coordinate (i.e. space) dimension.
idCounter which should be incremented to supply consistent vertex IDs.
justConfigIf true, then the configuration Element::m_conf will be updated but no nodes will be generated. This is used when considering boundary elements, which just require copying of face or edge interior nodes.

Reimplemented from Nektar::NekMeshUtils::Element.

Definition at line 285 of file Prism.cpp.

References ASSERTL1, Nektar::LibUtilities::PointsKey::GetPointsDim(), Nektar::NekMeshUtils::Element::m_conf, Nektar::NekMeshUtils::Element::m_curveType, Nektar::NekMeshUtils::ElmtConfig::m_faceNodes, Nektar::NekMeshUtils::ElmtConfig::m_order, Nektar::NekMeshUtils::ElmtConfig::m_volumeNodes, Nektar::NekMeshUtils::Element::m_volumeNodes, class_topology::Node, and Nektar::LibUtilities::PointsManager().

Referenced by ~Prism().

291 {
292  m_conf.m_order = order;
293  m_curveType = pType;
294  m_volumeNodes.clear();
295 
296  if (order == 1)
297  {
299  return;
300  }
301  else if (order == 2)
302  {
303  m_conf.m_faceNodes = true;
304  m_conf.m_volumeNodes = false;
305  return;
306  }
307 
308  m_conf.m_faceNodes = true;
309  m_conf.m_volumeNodes = true;
310 
311  if (justConfig)
312  {
313  return;
314  }
315 
316  int nPoints = order + 1;
317  StdRegions::StdExpansionSharedPtr xmap = geom->GetXmap();
318 
319  Array<OneD, NekDouble> px, py, pz;
320  LibUtilities::PointsKey pKey(nPoints, pType);
321  ASSERTL1(pKey.GetPointsDim() == 3, "Points distribution must be 3D");
322  LibUtilities::PointsManager()[pKey]->GetPoints(px, py, pz);
323 
324  Array<OneD, Array<OneD, NekDouble> > phys(coordDim);
325 
326  for (int i = 0; i < coordDim; ++i)
327  {
328  phys[i] = Array<OneD, NekDouble>(xmap->GetTotPoints());
329  xmap->BwdTrans(geom->GetCoeffs(i), phys[i]);
330  }
331 
332  const int nPrismPts = nPoints * nPoints * (nPoints + 1) / 2;
333  const int nPrismIntPts = (nPoints - 2) * (nPoints - 3) * (nPoints - 2) / 2;
334  m_volumeNodes.resize(nPrismIntPts);
335 
336  for (int i = nPrismPts - nPrismIntPts, cnt = 0; i < nPrismPts; ++i, ++cnt)
337  {
338  Array<OneD, NekDouble> xp(3);
339  xp[0] = px[i];
340  xp[1] = py[i];
341  xp[2] = pz[i];
342 
343  Array<OneD, NekDouble> x(3, 0.0);
344  for (int j = 0; j < coordDim; ++j)
345  {
346  x[j] = xmap->PhysEvaluate(xp, phys[j]);
347  }
348 
349  m_volumeNodes[cnt] = std::shared_ptr<Node>(
350  new Node(id++, x[0], x[1], x[2]));
351  }
352 }
bool m_faceNodes
Denotes whether the element contains face nodes. For 2D elements, if this is true then the element co...
Definition: ElementConfig.h:81
ElmtConfig m_conf
Contains configuration of the element.
Definition: Element.h:462
unsigned int m_order
Order of the element.
Definition: ElementConfig.h:88
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
bool m_volumeNodes
Denotes whether the element contains volume (i.e. interior) nodes. These are not supported by either ...
Definition: ElementConfig.h:86
PointsManagerT & PointsManager(void)
std::vector< NodeSharedPtr > m_volumeNodes
List of element volume nodes.
Definition: Element.h:474
LibUtilities::PointsType m_curveType
Volume curve type.
Definition: Element.h:476
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250

◆ OrientPrism()

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 523 of file Prism.cpp.

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

Referenced by Prism().

524 {
525  int lid[6], gid[6];
526 
527  // Re-order vertices.
528  for (int i = 0; i < 6; ++i)
529  {
530  lid[i] = i;
531  gid[i] = m_vertex[i]->m_id;
532  }
533 
534  gid[0] = gid[3] = max(gid[0], gid[3]);
535  gid[1] = gid[2] = max(gid[1], gid[2]);
536  gid[4] = gid[5] = max(gid[4], gid[5]);
537 
538  for (int i = 1; i < 6; ++i)
539  {
540  if (gid[0] < gid[i])
541  {
542  swap(gid[i], gid[0]);
543  swap(lid[i], lid[0]);
544  }
545  }
546 
547  if (lid[0] == 4 || lid[0] == 5)
548  {
549  m_orientation = 0;
550  }
551  else if (lid[0] == 1 || lid[0] == 2)
552  {
553  // Rotate prism clockwise in p-r plane
554  vector<NodeSharedPtr> vertexmap(6);
555  vertexmap[0] = m_vertex[4];
556  vertexmap[1] = m_vertex[0];
557  vertexmap[2] = m_vertex[3];
558  vertexmap[3] = m_vertex[5];
559  vertexmap[4] = m_vertex[1];
560  vertexmap[5] = m_vertex[2];
561  m_vertex = vertexmap;
562  m_orientation = 1;
563  }
564  else if (lid[0] == 0 || lid[0] == 3)
565  {
566  // Rotate prism counter-clockwise in p-r plane
567  vector<NodeSharedPtr> vertexmap(6);
568  vertexmap[0] = m_vertex[1];
569  vertexmap[1] = m_vertex[4];
570  vertexmap[2] = m_vertex[5];
571  vertexmap[3] = m_vertex[2];
572  vertexmap[4] = m_vertex[0];
573  vertexmap[5] = m_vertex[3];
574  m_vertex = vertexmap;
575  m_orientation = 2;
576  }
577  else
578  {
579  cerr << "Warning: possible prism orientation problem." << endl;
580  }
581 }
unsigned int m_orientation
Definition: Prism.h:94
std::vector< NodeSharedPtr > m_vertex
List of element vertex nodes.
Definition: Element.h:468

Member Data Documentation

◆ m_faceIds

int Nektar::NekMeshUtils::Prism::m_faceIds
staticprivate
Initial value:
= {
{0, 1, 2, 3}, {0, 1, 4, -1}, {1, 2, 5, 4}, {3, 2, 5, -1}, {0, 3, 5, 4}
}

Vertex IDs that make up prism faces.

Definition at line 100 of file Prism.h.

Referenced by GetFaceVertex(), and Prism().

◆ m_orientation

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 94 of file Prism.h.

Referenced by OrientPrism(), and Prism().

◆ m_type

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

Element type.

Definition at line 60 of file Prism.h.