Nektar++
Public Member Functions | Static Public Member Functions | Public Attributes | Static Public Attributes | Protected Member Functions | List of all members
Nektar::Utilities::Tetrahedron Class Reference

A 3-dimensional four-faced element. More...

#include <MeshElements.h>

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

Public Member Functions

 Tetrahedron (ElmtConfig pConf, std::vector< NodeSharedPtr > pNodeList, std::vector< int > pTagList)
 Create a tetrahedron element. More...
 
 Tetrahedron (const Tetrahedron &pSrc)
 
virtual ~Tetrahedron ()
 
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 tetrahedron. More...
 

Public Attributes

int m_orientationMap [4]
 
int m_origVertMap [4]
 

Static Public Attributes

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

Protected Member Functions

void OrientTet ()
 Orient tetrahedron 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 four-faced element.

Definition at line 1418 of file MeshElements.h.

Constructor & Destructor Documentation

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

Create a tetrahedron element.

Definition at line 840 of file MeshElements.cpp.

References Nektar::Utilities::HOTriangle< T >::Align(), ASSERTL0, 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_orientationMap, Nektar::Utilities::ElmtConfig::m_reorient, Nektar::Utilities::Element::m_tag, Nektar::Utilities::Element::m_taglist, Nektar::Utilities::Element::m_vertex, OrientTet(), and Nektar::Utilities::HOTriangle< T >::surfVerts.

843  : Element(pConf, GetNumNodes(pConf), pNodeList.size())
844  {
845  m_tag = "A";
846  m_dim = 3;
847  m_taglist = pTagList;
848  int n = m_conf.m_order-1;
849 
850  // Create a map to relate edge nodes to a pair of vertices
851  // defining an edge.
852  map<pair<int,int>, int> edgeNodeMap;
853  map<pair<int,int>, int>::iterator it;
854  edgeNodeMap[pair<int,int>(1,2)] = 5;
855  edgeNodeMap[pair<int,int>(2,3)] = 5 + n;
856  edgeNodeMap[pair<int,int>(1,3)] = 5 + 2*n;
857  edgeNodeMap[pair<int,int>(1,4)] = 5 + 3*n;
858  edgeNodeMap[pair<int,int>(2,4)] = 5 + 4*n;
859  edgeNodeMap[pair<int,int>(3,4)] = 5 + 5*n;
860 
861  // Add vertices
862  for (int i = 0; i < 4; ++i)
863  {
864  m_vertex.push_back(pNodeList[i]);
865  }
866 
867  // Create edges (with corresponding set of edge points)
868  int eid = 0;
869  for (it = edgeNodeMap.begin(); it != edgeNodeMap.end(); ++it)
870  {
871  vector<NodeSharedPtr> edgeNodes;
872  if (m_conf.m_order > 1) {
873  for (int j = it->second; j < it->second + n; ++j) {
874  edgeNodes.push_back(pNodeList[j-1]);
875  }
876  }
877  m_edge.push_back(
878  EdgeSharedPtr(new Edge(pNodeList[it->first.first-1],
879  pNodeList[it->first.second-1],
880  edgeNodes,
882  m_edge.back()->m_id = eid++;
883  }
884 
885  // Reorient the tet to ensure collapsed coordinates align between
886  // adjacent elements.
887  if (m_conf.m_reorient)
888  {
889  OrientTet();
890  }
891  else
892  {
893  // If we didn't need to orient the tet then set up the
894  // orientation map as the identity mapping.
895  for (int i = 0; i < 4; ++i)
896  {
897  m_orientationMap[i] = i;
898  }
899  }
900 
901  // Face-vertex IDs
902  int face_ids[4][3] = {
903  {0,1,2}, {0,1,3}, {1,2,3}, {0,2,3}
904  };
905 
906  // Face-edge IDs
907  int face_edges[4][3];
908 
909  // Create faces
910  for (int j = 0; j < 4; ++j)
911  {
912  vector<NodeSharedPtr> faceVertices;
913  vector<EdgeSharedPtr> faceEdges;
914  vector<NodeSharedPtr> faceNodes;
915 
916  // Extract the edges for this face. We need to do this because
917  // of the reorientation which might have been applied (see the
918  // additional note below).
919  for (int k = 0; k < 3; ++k)
920  {
921  faceVertices.push_back(m_vertex[face_ids[j][k]]);
922  NodeSharedPtr a = m_vertex[face_ids[j][k]];
923  NodeSharedPtr b = m_vertex[face_ids[j][(k+1)%3]];
924  for (unsigned int i = 0; i < m_edge.size(); ++i)
925  {
926  if ( ((*(m_edge[i]->m_n1)==*a) && (*(m_edge[i]->m_n2)==*b))
927  || ((*(m_edge[i]->m_n1)==*b) && (*(m_edge[i]->m_n2) == *a)) )
928  {
929  face_edges[j][k] = i;
930  faceEdges.push_back(m_edge[i]);
931  break;
932  }
933  }
934  }
935 
936  // When face curvature is supplied, it may have been the case
937  // that our tetrahedron was reoriented. In this case, faces have
938  // different vertex IDs and so we have to rotate the face
939  // curvature so that the two align appropriately.
940  if (m_conf.m_faceNodes)
941  {
942  const int nFaceNodes = n*(n-1)/2;
943 
944  // Get the vertex IDs of whatever face we are processing.
945  vector<int> faceIds(3);
946  faceIds[0] = faceVertices[0]->m_id;
947  faceIds[1] = faceVertices[1]->m_id;
948  faceIds[2] = faceVertices[2]->m_id;
949 
950  // Find out the original face number as we were given it in
951  // the constructor using the orientation map.
952  int origFace = -1;
953  for (int i = 0; i < 4; ++i)
954  {
955  if (m_orientationMap[i] == j)
956  {
957  origFace = i;
958  break;
959  }
960  }
961 
962  ASSERTL0(origFace >= 0, "Couldn't find face");
963 
964  // Now get the face nodes for the original face.
965  int N = 4 + 6*n + origFace * nFaceNodes;
966  for (int i = 0; i < nFaceNodes; ++i)
967  {
968  faceNodes.push_back(pNodeList[N+i]);
969  }
970 
971  // Find the original face vertex IDs.
972  vector<int> origFaceIds(3);
973  origFaceIds[0] = pNodeList[face_ids[origFace][0]]->m_id;
974  origFaceIds[1] = pNodeList[face_ids[origFace][1]]->m_id;
975  origFaceIds[2] = pNodeList[face_ids[origFace][2]]->m_id;
976 
977  // Construct a HOTriangle object which performs the
978  // orientation magically for us.
979  HOTriangle<NodeSharedPtr> hoTri(origFaceIds, faceNodes);
980  hoTri.Align(faceIds);
981 
982  // Copy the face nodes back again.
983  faceNodes = hoTri.surfVerts;
984  }
985 
986  m_face.push_back(FaceSharedPtr(
987  new Face(faceVertices, faceNodes, faceEdges, m_conf.m_faceCurveType)));
988  }
989 
990  vector<EdgeSharedPtr> tmp(6);
991  tmp[0] = m_edge[face_edges[0][0]];
992  tmp[1] = m_edge[face_edges[0][1]];
993  tmp[2] = m_edge[face_edges[0][2]];
994  tmp[3] = m_edge[face_edges[1][2]];
995  tmp[4] = m_edge[face_edges[1][1]];
996  tmp[5] = m_edge[face_edges[2][1]];
997  m_edge = tmp;
998  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
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
static unsigned int GetNumNodes(ElmtConfig pConf)
Return the number of nodes defining a tetrahedron.
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
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
void OrientTet()
Orient tetrahedron to align degenerate vertices.
Element(ElmtConfig pConf, unsigned int pNumNodes, unsigned int pGotNodes)
ElmtConfig m_conf
Contains configuration of the element.
Definition: MeshElements.h:980
Nektar::Utilities::Tetrahedron::Tetrahedron ( const Tetrahedron pSrc)
virtual Nektar::Utilities::Tetrahedron::~Tetrahedron ( )
inlinevirtual

Definition at line 1442 of file MeshElements.h.

1442 {}

Member Function Documentation

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

Reimplemented from Nektar::Utilities::Element.

Definition at line 1034 of file MeshElements.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eGaussRadauMAlpha2Beta0, Nektar::LibUtilities::eNodalTetEvenlySpaced, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eOrtho_C, 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.

1035  {
1036  int i, j;
1037 
1038  // Create basis key for a nodal tetrahedron.
1039  LibUtilities::BasisKey B0(
1040  LibUtilities::eOrtho_A, order+1,
1041  LibUtilities::PointsKey(
1043  LibUtilities::BasisKey B1(
1044  LibUtilities::eOrtho_B, order+1,
1045  LibUtilities::PointsKey(
1047  LibUtilities::BasisKey B2(
1048  LibUtilities::eOrtho_C, order+1,
1049  LibUtilities::PointsKey(
1051 
1052  // Create a standard nodal tetrahedron in order to get the
1053  // Vandermonde matrix to perform interpolation to nodal points.
1057 
1058  Array<OneD, NekDouble> x, y, z;
1059 
1060  nodalTet->GetNodalPoints(x,y,z);
1061 
1063  boost::dynamic_pointer_cast<SpatialDomains::TetGeom>(
1064  this->GetGeom(3));
1065 
1066  // Create basis key for a tetrahedron.
1067  LibUtilities::BasisKey C0(
1068  LibUtilities::eOrtho_A, order+1,
1069  LibUtilities::PointsKey(
1071  LibUtilities::BasisKey C1(
1072  LibUtilities::eOrtho_B, order+1,
1073  LibUtilities::PointsKey(
1075  LibUtilities::BasisKey C2(
1076  LibUtilities::eOrtho_C, order+1,
1077  LibUtilities::PointsKey(
1079 
1080  // Create a tet.
1083  C0, C1, C2, geom);
1084 
1085  // Get coordinate array for tetrahedron.
1086  int nqtot = tet->GetTotPoints();
1087  Array<OneD, NekDouble> alloc(6*nqtot);
1088  Array<OneD, NekDouble> xi(alloc);
1089  Array<OneD, NekDouble> yi(alloc+ nqtot);
1090  Array<OneD, NekDouble> zi(alloc+2*nqtot);
1091  Array<OneD, NekDouble> xo(alloc+3*nqtot);
1092  Array<OneD, NekDouble> yo(alloc+4*nqtot);
1093  Array<OneD, NekDouble> zo(alloc+5*nqtot);
1094  Array<OneD, NekDouble> tmp;
1095 
1096  tet->GetCoords(xi, yi, zi);
1097 
1098  for (i = 0; i < 3; ++i)
1099  {
1100  Array<OneD, NekDouble> coeffs(nodalTet->GetNcoeffs());
1101  tet->FwdTrans(alloc+i*nqtot, coeffs);
1102  // Apply Vandermonde matrix to project onto nodal space.
1103  nodalTet->ModalToNodal(coeffs, tmp=alloc+(i+3)*nqtot);
1104  }
1105 
1106  // Now extract points from the co-ordinate arrays into the
1107  // edge/face/volume nodes. First, extract edge-interior nodes.
1108  for (i = 0; i < 6; ++i)
1109  {
1110  int pos = 4 + i*(order-1);
1111  m_edge[i]->m_edgeNodes.clear();
1112  for (j = 0; j < order-1; ++j)
1113  {
1114  m_edge[i]->m_edgeNodes.push_back(
1115  NodeSharedPtr(new Node(0, xo[pos+j], yo[pos+j], zo[pos+j])));
1116  }
1117  }
1118 
1119  // Now extract face-interior nodes.
1120  for (i = 0; i < 4; ++i)
1121  {
1122  int pos = 4 + 6*(order-1) + i*(order-2)*(order-1)/2;
1123  m_face[i]->m_faceNodes.clear();
1124  for (j = 0; j < (order-2)*(order-1)/2; ++j)
1125  {
1126  m_face[i]->m_faceNodes.push_back(
1127  NodeSharedPtr(new Node(0, xo[pos+j], yo[pos+j], zo[pos+j])));
1128  }
1129  }
1130 
1131  // Finally extract volume nodes.
1132  int pos = 4 + 6*(order-1) + 4*(order-2)*(order-1)/2;
1133  for (i = pos; i < (order+1)*(order+2)*(order+3)/6; ++i)
1134  {
1135  m_volumeNodes.push_back(
1136  NodeSharedPtr(new Node(0, xo[i], yo[i], zo[i])));
1137  }
1138 
1139  m_conf.m_order = order;
1140  m_conf.m_faceNodes = true;
1141  m_conf.m_volumeNodes = true;
1142  }
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
virtual SpatialDomains::GeometrySharedPtr GetGeom(int coordDim)
Generate a Nektar++ geometry object for this element.
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
boost::shared_ptr< TetExp > TetExpSharedPtr
Definition: TetExp.h:234
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:57
Principle Orthogonal Functions .
Definition: BasisType.h:47
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:48
Principle Orthogonal Functions .
Definition: BasisType.h:46
boost::shared_ptr< StdNodalTetExp > StdNodalTetExpSharedPtr
3D Evenly-spaced points on a Tetrahedron
Definition: PointsType.h:71
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
boost::shared_ptr< TetGeom > TetGeomSharedPtr
Definition: TetGeom.h:106
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::Tetrahedron::create ( ElmtConfig  pConf,
std::vector< NodeSharedPtr pNodeList,
std::vector< int >  pTagList 
)
inlinestatic

Creates an instance of this class.

Definition at line 1421 of file MeshElements.h.

1425  {
1426  ElementSharedPtr e = boost::shared_ptr<Element>(
1427  new Tetrahedron(pConf, pNodeList, pTagList));
1428  vector<FaceSharedPtr> faces = e->GetFaceList();
1429  for (int i = 0; i < faces.size(); ++i)
1430  {
1431  faces[i]->m_elLink.push_back(pair<ElementSharedPtr, int>(e,i));
1432  }
1433  return e;
1434  }
boost::shared_ptr< Element > ElementSharedPtr
Shared pointer to an element.
Definition: MeshElements.h:63
Tetrahedron(ElmtConfig pConf, std::vector< NodeSharedPtr > pNodeList, std::vector< int > pTagList)
Create a tetrahedron element.
SpatialDomains::GeometrySharedPtr Nektar::Utilities::Tetrahedron::GetGeom ( int  coordDim)
virtual

Generate a Nektar++ geometry object for this element.

Reimplemented from Nektar::Utilities::Element.

Definition at line 1000 of file MeshElements.cpp.

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

Referenced by Complete().

1001  {
1004 
1005  for (int i = 0; i < 4; ++i)
1006  {
1007  tfaces[i] = boost::dynamic_pointer_cast
1008  <SpatialDomains::TriGeom>(m_face[i]->GetGeom(coordDim));
1009  }
1010 
1012  AllocateSharedPtr(tfaces);
1013 
1014  return ret;
1015  }
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< TetGeom > TetGeomSharedPtr
Definition: TetGeom.h:106
boost::shared_ptr< TriGeom > TriGeomSharedPtr
Definition: TriGeom.h:58
unsigned int Nektar::Utilities::Tetrahedron::GetNumNodes ( ElmtConfig  pConf)
static

Return the number of nodes defining a tetrahedron.

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

1021  {
1022  int n = pConf.m_order;
1023  if (pConf.m_volumeNodes && pConf.m_faceNodes)
1024  return (n+1)*(n+2)*(n+3)/6;
1025  else if (!pConf.m_volumeNodes && pConf.m_faceNodes)
1026  return 4*(n+1)*(n+2)/2-6*(n+1)+4;
1027  else
1028  return 6*(n+1)-8;
1029  }
void Nektar::Utilities::Tetrahedron::OrientTet ( )
protected

Orient tetrahedron to align degenerate vertices.

Orientation of tetrahedral 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 vertices are ordered with the highest global vertex at the top degenerate point, and the base degenerate point has second lowest ID. These vertices are swapped if the element is incorrectly oriented.

Definition at line 1191 of file MeshElements.cpp.

References Nektar::iterator, Nektar::Utilities::Element::m_id, m_orientationMap, m_origVertMap, and Nektar::Utilities::Element::m_vertex.

Referenced by Tetrahedron().

1192  {
1193  static int face_ids[4][3] = {
1194  {0,1,2},{0,1,3},{1,2,3},{0,2,3}};
1195 
1196  TetOrientSet faces;
1197 
1198  // Create a copy of the original vertex ordering. This is used to
1199  // construct a mapping, #orientationMap, which maps the original
1200  // face ordering to the new face ordering.
1201  for (int i = 0; i < 4; ++i)
1202  {
1203  vector<int> nodes(3);
1204 
1205  nodes[0] = m_vertex[face_ids[i][0]]->m_id;
1206  nodes[1] = m_vertex[face_ids[i][1]]->m_id;
1207  nodes[2] = m_vertex[face_ids[i][2]]->m_id;
1208 
1209  sort(nodes.begin(), nodes.end());
1210  struct TetOrient faceNodes(nodes, i);
1211  faces.insert(faceNodes);
1212  }
1213 
1214  // Store a copy of the original vertex ordering so we can create a
1215  // permutation map later.
1216  vector<NodeSharedPtr> origVert = m_vertex;
1217 
1218  // Order vertices with highest global vertex at top degenerate
1219  // point. Place second highest global vertex at base degenerate
1220  // point.
1221  sort(m_vertex.begin(), m_vertex.end());
1222 
1223  // Calculate a.(b x c) if negative, reverse order of
1224  // non-degenerate points to correctly orientate the tet.
1225 
1226  NekDouble ax = m_vertex[1]->m_x-m_vertex[0]->m_x;
1227  NekDouble ay = m_vertex[1]->m_y-m_vertex[0]->m_y;
1228  NekDouble az = m_vertex[1]->m_z-m_vertex[0]->m_z;
1229  NekDouble bx = m_vertex[2]->m_x-m_vertex[0]->m_x;
1230  NekDouble by = m_vertex[2]->m_y-m_vertex[0]->m_y;
1231  NekDouble bz = m_vertex[2]->m_z-m_vertex[0]->m_z;
1232  NekDouble cx = m_vertex[3]->m_x-m_vertex[0]->m_x;
1233  NekDouble cy = m_vertex[3]->m_y-m_vertex[0]->m_y;
1234  NekDouble cz = m_vertex[3]->m_z-m_vertex[0]->m_z;
1235 
1236  NekDouble nx = (ay*bz-az*by);
1237  NekDouble ny = (az*bx-ax*bz);
1238  NekDouble nz = (ax*by-ay*bx);
1239  NekDouble nmag = sqrt(nx*nx+ny*ny+nz*nz);
1240  nx /= nmag; ny /= nmag; nz /= nmag;
1241 
1242  // distance of top vertex from base
1243  NekDouble dist = cx*nx+cy*ny+cz*nz;
1244 
1245  if (fabs(dist) <= 1e-4)
1246  {
1247  cerr << "Warning: degenerate tetrahedron, 3rd vertex is = " << dist <<" from face" << endl;
1248  }
1249 
1250  if (dist < 0)
1251  {
1252  swap(m_vertex[0], m_vertex[1]);
1253  }
1254 
1255  nx = (ay*cz-az*cy);
1256  ny = (az*cx-ax*cz);
1257  nz = (ax*cy-ay*cx);
1258  nmag = sqrt(nx*nx+ny*ny+nz*nz);
1259  nx /= nmag; ny /= nmag; nz /= nmag;
1260 
1261  // distance of top vertex from base
1262  dist = bx*nx+by*ny+bz*nz;
1263 
1264  if (fabs(dist) <= 1e-4)
1265  {
1266  cerr << "Warning: degenerate tetrahedron, 2nd vertex is = " << dist <<" from face" << endl;
1267  }
1268 
1269  nx = (by*cz-bz*cy);
1270  ny = (bz*cx-bx*cz);
1271  nz = (bx*cy-by*cx);
1272  nmag = sqrt(nx*nx+ny*ny+nz*nz);
1273  nx /= nmag; ny /= nmag; nz /= nmag;
1274 
1275  // distance of top vertex from base
1276  dist = ax*nx+ay*ny+az*nz;
1277 
1278  if (fabs(dist) <= 1e-4)
1279  {
1280  cerr << "Warning: degenerate tetrahedron, 1st vertex is = " << dist <<" from face" << endl;
1281  }
1282 
1283 
1285 
1286  // Search for the face in the original set of face nodes. Then use
1287  // this to construct the #orientationMap.
1288  for (int i = 0; i < 4; ++i)
1289  {
1290  vector<int> nodes(3);
1291 
1292  nodes[0] = m_vertex[face_ids[i][0]]->m_id;
1293  nodes[1] = m_vertex[face_ids[i][1]]->m_id;
1294  nodes[2] = m_vertex[face_ids[i][2]]->m_id;
1295  sort(nodes.begin(), nodes.end());
1296 
1297  struct TetOrient faceNodes(nodes, 0);
1298 
1299  it = faces.find(faceNodes);
1300  m_orientationMap[it->fid] = i;
1301 
1302  for (int j = 0; j < 4; ++j)
1303  {
1304  if (m_vertex[i]->m_id == origVert[j]->m_id)
1305  {
1306  m_origVertMap[j] = i;
1307  break;
1308  }
1309  }
1310  }
1311  }
unsigned int m_id
ID of the element.
Definition: MeshElements.h:976
double NekDouble
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
std::vector< NodeSharedPtr > m_vertex
List of element vertex nodes.
Definition: MeshElements.h:986
boost::unordered_set< struct TetOrient, TetOrientHash > TetOrientSet

Member Data Documentation

int Nektar::Utilities::Tetrahedron::m_orientationMap[4]

Definition at line 1449 of file MeshElements.h.

Referenced by OrientTet(), Nektar::Utilities::InputNek::Process(), and Tetrahedron().

int Nektar::Utilities::Tetrahedron::m_origVertMap[4]

Definition at line 1450 of file MeshElements.h.

Referenced by OrientTet().

LibUtilities::ShapeType Nektar::Utilities::Tetrahedron::m_type
static
Initial value:
RegisterCreatorFunction(LibUtilities::eTetrahedron, Tetrahedron::create, "Tetrahedron")

Element type.

Definition at line 1436 of file MeshElements.h.