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

Static Public Member Functions

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

Public Attributes

int m_orientationMap [4]

Static Public Attributes

static LibUtilities::ShapeType m_type
 Element type.

Protected Member Functions

void OrientTet ()
 Orient tetrahedron to align degenerate vertices.

Additional Inherited Members

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

Detailed Description

A 3-dimensional four-faced element.

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

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

: Element(pConf, GetNumNodes(pConf), pNodeList.size())
{
m_tag = "A";
m_dim = 3;
m_taglist = pTagList;
int n = m_conf.m_order-1;
// Create a map to relate edge nodes to a pair of vertices
// defining an edge.
map<pair<int,int>, int> edgeNodeMap;
map<pair<int,int>, int>::iterator it;
edgeNodeMap[pair<int,int>(1,2)] = 5;
edgeNodeMap[pair<int,int>(2,3)] = 5 + n;
edgeNodeMap[pair<int,int>(1,3)] = 5 + 2*n;
edgeNodeMap[pair<int,int>(1,4)] = 5 + 3*n;
edgeNodeMap[pair<int,int>(2,4)] = 5 + 4*n;
edgeNodeMap[pair<int,int>(3,4)] = 5 + 5*n;
// Add vertices
for (int i = 0; i < 4; ++i) {
m_vertex.push_back(pNodeList[i]);
}
// Create edges (with corresponding set of edge points)
int eid = 0;
for (it = edgeNodeMap.begin(); it != edgeNodeMap.end(); ++it)
{
vector<NodeSharedPtr> edgeNodes;
if (m_conf.m_order > 1) {
for (int j = it->second; j < it->second + n; ++j) {
edgeNodes.push_back(pNodeList[j-1]);
}
}
m_edge.push_back(EdgeSharedPtr(new Edge(pNodeList[it->first.first-1],
pNodeList[it->first.second-1],
edgeNodes,
m_edge.back()->m_id = eid++;
}
// Reorient the tet to ensure collapsed coordinates align between adjacent
// elements.
{
}
// Create faces
int face_ids[4][3] = {
{0,1,2},{0,1,3},{1,2,3},{0,2,3}};
int face_edges[4][3];
for (int j = 0; j < 4; ++j)
{
vector<NodeSharedPtr> faceVertices;
vector<EdgeSharedPtr> faceEdges;
vector<NodeSharedPtr> faceNodes;
for (int k = 0; k < 3; ++k)
{
faceVertices.push_back(m_vertex[face_ids[j][k]]);
NodeSharedPtr a = m_vertex[face_ids[j][k]];
NodeSharedPtr b = m_vertex[face_ids[j][(k+1)%3]];
for (unsigned int i = 0; i < m_edge.size(); ++i)
{
if ( ((*(m_edge[i]->m_n1)==*a) && (*(m_edge[i]->m_n2)==*b))
|| ((*(m_edge[i]->m_n1)==*b) && (*(m_edge[i]->m_n2) == *a)) )
{
face_edges[j][k] = i;
faceEdges.push_back(m_edge[i]);
break;
}
}
}
{
int N = 4 + 6*n + j*n*(n-1)/2;
for (int i = 0; i < n*(n-1)/2; ++i)
{
faceNodes.push_back(pNodeList[N+i]);
}
}
m_face.push_back(FaceSharedPtr(
new Face(faceVertices, faceNodes, faceEdges, m_conf.m_faceCurveType)));
}
vector<EdgeSharedPtr> tmp(6);
tmp[0] = m_edge[face_edges[0][0]];
tmp[1] = m_edge[face_edges[0][1]];
tmp[2] = m_edge[face_edges[0][2]];
tmp[3] = m_edge[face_edges[1][2]];
tmp[4] = m_edge[face_edges[1][1]];
tmp[5] = m_edge[face_edges[2][1]];
m_edge = tmp;
}
Nektar::Utilities::Tetrahedron::Tetrahedron ( const Tetrahedron pSrc)
virtual Nektar::Utilities::Tetrahedron::~Tetrahedron ( )
inlinevirtual

Definition at line 1363 of file MeshElements.h.

{}

Member Function Documentation

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

Reimplemented from Nektar::Utilities::Element.

Definition at line 980 of file MeshElements.cpp.

References 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.

{
int i, j;
// Create basis key for a nodal tetrahedron.
LibUtilities::BasisKey B0(
LibUtilities::PointsKey(
LibUtilities::BasisKey B1(
LibUtilities::PointsKey(
LibUtilities::BasisKey B2(
LibUtilities::PointsKey(
// Create a standard nodal tetrahedron in order to get the
// Vandermonde matrix to perform interpolation to nodal points.
MemoryManager<StdRegions::StdNodalTetExp>::AllocateSharedPtr(
Array<OneD, NekDouble> x, y, z;
nodalTet->GetNodalPoints(x,y,z);
boost::dynamic_pointer_cast<SpatialDomains::TetGeom>(
this->GetGeom(3));
// Create basis key for a tetrahedron.
LibUtilities::BasisKey C0(
LibUtilities::PointsKey(
LibUtilities::BasisKey C1(
LibUtilities::PointsKey(
LibUtilities::BasisKey C2(
LibUtilities::PointsKey(
// Create a tet.
MemoryManager<LocalRegions::TetExp>::AllocateSharedPtr(
C0, C1, C2, geom);
// Get coordinate array for tetrahedron.
int nqtot = tet->GetTotPoints();
Array<OneD, NekDouble> alloc(6*nqtot);
Array<OneD, NekDouble> xi(alloc);
Array<OneD, NekDouble> yi(alloc+ nqtot);
Array<OneD, NekDouble> zi(alloc+2*nqtot);
Array<OneD, NekDouble> xo(alloc+3*nqtot);
Array<OneD, NekDouble> yo(alloc+4*nqtot);
Array<OneD, NekDouble> zo(alloc+5*nqtot);
Array<OneD, NekDouble> tmp;
tet->GetCoords(xi, yi, zi);
for (i = 0; i < 3; ++i)
{
Array<OneD, NekDouble> coeffs(nodalTet->GetNcoeffs());
tet->FwdTrans(alloc+i*nqtot, coeffs);
// Apply Vandermonde matrix to project onto nodal space.
nodalTet->ModalToNodal(coeffs, tmp=alloc+(i+3)*nqtot);
}
// Now extract points from the co-ordinate arrays into the
// edge/face/volume nodes. First, extract edge-interior nodes.
for (i = 0; i < 6; ++i)
{
int pos = 4 + i*(order-1);
m_edge[i]->m_edgeNodes.clear();
for (j = 0; j < order-1; ++j)
{
m_edge[i]->m_edgeNodes.push_back(
NodeSharedPtr(new Node(0, xo[pos+j], yo[pos+j], zo[pos+j])));
}
}
// Now extract face-interior nodes.
for (i = 0; i < 4; ++i)
{
int pos = 4 + 6*(order-1) + i*(order-2)*(order-1)/2;
m_face[i]->m_faceNodes.clear();
for (j = 0; j < (order-2)*(order-1)/2; ++j)
{
m_face[i]->m_faceNodes.push_back(
NodeSharedPtr(new Node(0, xo[pos+j], yo[pos+j], zo[pos+j])));
}
}
// Finally extract volume nodes.
int pos = 4 + 6*(order-1) + 4*(order-2)*(order-1)/2;
for (i = pos; i < (order+1)*(order+2)*(order+3)/6; ++i)
{
m_volumeNodes.push_back(
NodeSharedPtr(new Node(0, xo[i], yo[i], zo[i])));
}
m_conf.m_order = order;
}
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 1342 of file MeshElements.h.

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

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

Generate a Nektar++ geometry object for this element.

Reimplemented from Nektar::Utilities::Element.

Definition at line 946 of file MeshElements.cpp.

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

Referenced by Complete().

{
for (int i = 0; i < 4; ++i)
{
tfaces[i] = boost::dynamic_pointer_cast
<SpatialDomains::TriGeom>(m_face[i]->GetGeom(coordDim));
}
ret = MemoryManager<SpatialDomains::TetGeom>::
AllocateSharedPtr(tfaces);
return ret;
}
unsigned int Nektar::Utilities::Tetrahedron::GetNumNodes ( ElmtConfig  pConf)
static

Return the number of nodes defining a tetrahedron.

Definition at line 966 of file MeshElements.cpp.

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

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

{
int n = pConf.m_order;
if (pConf.m_volumeNodes && pConf.m_faceNodes)
return (n+1)*(n+2)*(n+3)/6;
else if (!pConf.m_volumeNodes && pConf.m_faceNodes)
return 4*(n+1)*(n+2)/2-6*(n+1)+4;
else
return 6*(n+1)-8;
}
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 1137 of file MeshElements.cpp.

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

Referenced by Tetrahedron().

{
static int face_ids[4][3] = {
{0,1,2},{0,1,3},{1,2,3},{0,2,3}};
TetOrientSet faces;
// Create a copy of the original vertex ordering. This is used to
// construct a mapping, #orientationMap, which maps the original
// face ordering to the new face ordering.
for (int i = 0; i < 4; ++i)
{
vector<int> nodes(3);
nodes[0] = m_vertex[face_ids[i][0]]->m_id;
nodes[1] = m_vertex[face_ids[i][1]]->m_id;
nodes[2] = m_vertex[face_ids[i][2]]->m_id;
sort(nodes.begin(), nodes.end());
struct TetOrient faceNodes(nodes, i);
faces.insert(faceNodes);
}
// Order vertices with highest global vertex at top degenerate
// point. Place second highest global vertex at base degenerate
// point.
sort(m_vertex.begin(), m_vertex.end());
// Calculate a.(b x c) to determine tet volume; if negative,
// reverse order of non-degenerate points to correctly orientate
// the tet.
double ax = m_vertex[1]->m_x-m_vertex[0]->m_x;
double ay = m_vertex[1]->m_y-m_vertex[0]->m_y;
double az = m_vertex[1]->m_z-m_vertex[0]->m_z;
double bx = m_vertex[2]->m_x-m_vertex[0]->m_x;
double by = m_vertex[2]->m_y-m_vertex[0]->m_y;
double bz = m_vertex[2]->m_z-m_vertex[0]->m_z;
double cx = m_vertex[3]->m_x-m_vertex[0]->m_x;
double cy = m_vertex[3]->m_y-m_vertex[0]->m_y;
double cz = m_vertex[3]->m_z-m_vertex[0]->m_z;
double vol = cx*(ay*bz-az*by)+cy*(az*bx-ax*bz)+cz*(ax*by-ay*bx);
vol /= 6.0;
if (fabs(vol) <= 1e-10)
{
cerr << "Warning: degenerate tetrahedron, volume = " << vol << endl;
}
if (vol < 0)
{
swap(m_vertex[0], m_vertex[1]);
}
// Search for the face in the original set of face nodes. Then use
// this to construct the #orientationMap.
for (int i = 0; i < 4; ++i)
{
vector<int> nodes(3);
nodes[0] = m_vertex[face_ids[i][0]]->m_id;
nodes[1] = m_vertex[face_ids[i][1]]->m_id;
nodes[2] = m_vertex[face_ids[i][2]]->m_id;
sort(nodes.begin(), nodes.end());
struct TetOrient faceNodes(nodes, 0);
it = faces.find(faceNodes);
m_orientationMap[it->fid] = i;
}
}

Member Data Documentation

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

Orientation of tet; unchanged = 0; base vertex swapped = 1.

Definition at line 1373 of file MeshElements.h.

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

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

Element type.

Definition at line 1357 of file MeshElements.h.