Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
Nektar::SpatialDomains::TriGeom Class Reference

#include <TriGeom.h>

Inheritance diagram for Nektar::SpatialDomains::TriGeom:
Inheritance graph
[legend]
Collaboration diagram for Nektar::SpatialDomains::TriGeom:
Collaboration graph
[legend]

Public Member Functions

 TriGeom ()
 TriGeom (int id, const int coordim)
 TriGeom (const int id, const PointGeomSharedPtr verts[], const SegGeomSharedPtr edges[], const StdRegions::Orientation eorient[])
 TriGeom (const int id, const SegGeomSharedPtr edges[], const StdRegions::Orientation eorient[])
 TriGeom (const int id, const SegGeomSharedPtr edges[], const StdRegions::Orientation eorient[], const CurveSharedPtr &curve)
 TriGeom (const TriGeom &in)
 ~TriGeom ()
NekDouble GetCoord (const int i, const Array< OneD, const NekDouble > &Lcoord)
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry2D
 Geometry2D ()
 Geometry2D (const int coordim)
virtual ~Geometry2D ()
int GetFid () const
const Geometry1DSharedPtr GetEdge (int i) const
const Geometry2DSharedPtr GetFace (int i) const
int WhichEdge (SegGeomSharedPtr edge)
int WhichFace (Geometry2DSharedPtr face)
const LibUtilities::BasisSharedPtr GetEdgeBasis (const int i)
StdRegions::Orientation GetFaceOrient (const int i) const
StdRegions::Orientation GetCartesianEorient (const int i) const
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry
 Geometry ()
 Geometry (int coordim)
virtual ~Geometry ()
bool IsElmtConnected (int gvo_id, int locid) const
void AddElmtConnected (int gvo_id, int locid)
int NumElmtConnected () const
int GetCoordim () const
void SetCoordim (int coordim)
GeomFactorsSharedPtr GetGeomFactors ()
GeomFactorsSharedPtr GetRefGeomFactors (const Array< OneD, const LibUtilities::BasisSharedPtr > &tbasis)
GeomFactorsSharedPtr GetMetricInfo ()
LibUtilities::ShapeType GetShapeType (void)
int GetGlobalID (void)
void SetGlobalID (int globalid)
int GetVid (int i) const
int GetEid (int i) const
int GetFid (int i) const
int GetTid (int i) const
int GetNumVerts () const
PointGeomSharedPtr GetVertex (int i) const
StdRegions::Orientation GetEorient (const int i) const
StdRegions::Orientation GetPorient (const int i) const
StdRegions::Orientation GetForient (const int i) const
int GetNumEdges () const
int GetNumFaces () const
int GetShapeDim () const
StdRegions::StdExpansionSharedPtr GetXmap () const
const Array< OneD, const
NekDouble > & 
GetCoeffs (const int i) const
bool ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
bool ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol)
bool ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &resid)
int GetVertexEdgeMap (int i, int j) const
int GetVertexFaceMap (int i, int j) const
 return the id of the $j^{th}$ face attached to the $ i^{th}$ vertex
int GetEdgeFaceMap (int i, int j) const
void FillGeom ()
 Put all quadrature information into face/edge structure and backward transform.
NekDouble GetLocCoords (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
void SetOwnData ()
const LibUtilities::BasisSharedPtr GetBasis (const int i)
 Return the j-th basis of the i-th co-ordinate dimension.
const LibUtilities::PointsKeyVector GetPointsKeys ()

Static Public Member Functions

static StdRegions::Orientation GetFaceOrientation (const TriGeom &face1, const TriGeom &face2)
static StdRegions::Orientation GetFaceOrientation (const PointGeomVector &face1, const PointGeomVector &face2)

Static Public Attributes

static const int kNedges = 3
 Get the orientation of face1.
static const int kNverts = 3

Protected Member Functions

virtual void v_AddElmtConnected (int gvo_id, int locid)
virtual int v_NumElmtConnected () const
virtual bool v_IsElmtConnected (int gvo_id, int locid) const
virtual int v_GetFid () const
virtual int v_GetCoordim () const
virtual const
LibUtilities::BasisSharedPtr 
v_GetBasis (const int i)
virtual const
LibUtilities::BasisSharedPtr 
v_GetEdgeBasis (const int i)
virtual NekDouble v_GetCoord (const int i, const Array< OneD, const NekDouble > &Lcoord)
virtual void v_GenGeomFactors ()
virtual void v_SetOwnData ()
virtual void v_FillGeom ()
 Put all quadrature information into edge structure.
virtual NekDouble v_GetLocCoords (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
virtual int v_GetEid (int i) const
virtual int v_GetVid (int i) const
virtual PointGeomSharedPtr v_GetVertex (int i) const
virtual const Geometry1DSharedPtr v_GetEdge (int i) const
virtual StdRegions::Orientation v_GetEorient (const int i) const
virtual StdRegions::Orientation v_GetCartesianEorient (const int i) const
virtual int v_WhichEdge (SegGeomSharedPtr edge)
 Return the edge number of the given edge.
virtual int v_GetNumVerts () const
virtual int v_GetNumEdges () const
virtual bool v_ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
 Determines if a point specified in global coordinates is located within this tetrahedral geometry.
virtual bool v_ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol)
virtual bool v_ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &resid)
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry2D
void NewtonIterationForLocCoord (const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &ptsx, const Array< OneD, const NekDouble > &ptsy, Array< OneD, NekDouble > &Lcoords, NekDouble &resid)
- Protected Member Functions inherited from Nektar::SpatialDomains::Geometry
void GenGeomFactors ()
virtual int v_GetFid (int i) const
virtual StdRegions::Orientation v_GetPorient (const int i) const
virtual StdRegions::Orientation v_GetForient (const int i) const
virtual int v_GetNumFaces () const
virtual int v_GetShapeDim () const
virtual
StdRegions::StdExpansionSharedPtr 
v_GetXmap () const
virtual int v_GetVertexEdgeMap (int i, int j) const
virtual int v_GetVertexFaceMap (int i, int j) const
virtual int v_GetEdgeFaceMap (int i, int j) const
void SetUpCoeffs (const int nCoeffs)

Protected Attributes

PointGeomVector m_verts
SegGeomVector m_edges
StdRegions::Orientation m_eorient [kNedges]
int m_fid
bool m_ownVerts
std::list< CompToElmtm_elmtMap

Private Attributes

bool m_ownData

Additional Inherited Members

- Static Protected Member Functions inherited from Nektar::SpatialDomains::Geometry
static GeomFactorsSharedPtr ValidateRegGeomFactor (GeomFactorsSharedPtr geomFactor)
- Static Protected Attributes inherited from Nektar::SpatialDomains::Geometry
static GeomFactorsVector m_regGeomFactorsManager

Detailed Description

Definition at line 65 of file TriGeom.h.

Constructor & Destructor Documentation

Nektar::SpatialDomains::TriGeom::TriGeom ( )
Nektar::SpatialDomains::TriGeom::TriGeom ( int  id,
const int  coordim 
)
Nektar::SpatialDomains::TriGeom::TriGeom ( const int  id,
const PointGeomSharedPtr  verts[],
const SegGeomSharedPtr  edges[],
const StdRegions::Orientation  eorient[] 
)

Copy the vert shared pointers.

Copy the edge shared pointers.

Definition at line 79 of file TriGeom.cpp.

References ASSERTL0, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eTriangle, Nektar::SpatialDomains::Geometry::GetBasis(), kNedges, kNverts, Nektar::SpatialDomains::Geometry::m_coordim, m_edges, m_eorient, m_fid, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_shapeType, m_verts, Nektar::SpatialDomains::Geometry::m_xmap, and Nektar::SpatialDomains::Geometry::SetUpCoeffs().

:
Geometry2D(verts[0]->GetCoordim()),
m_fid(id)
{
/// Copy the vert shared pointers.
m_verts.insert(m_verts.begin(), verts, verts+TriGeom::kNverts);
/// Copy the edge shared pointers.
m_edges.insert(m_edges.begin(), edges, edges+TriGeom::kNedges);
for (int j=0; j<kNedges; ++j)
{
m_eorient[j] = eorient[j];
}
m_coordim = verts[0]->GetCoordim();
"Cannot call function with dim == 1");
int order0 = edges[0]->GetBasis(0)->GetNumModes();
int points0 = edges[0]->GetBasis(0)->GetNumPoints();
int order1 = max(order0,
max(edges[1]->GetBasis(0)->GetNumModes(),
edges[2]->GetBasis(0)->GetNumModes()));
int points1 = max(points0,
max(edges[1]->GetBasis(0)->GetNumPoints(),
edges[2]->GetBasis(0)->GetNumPoints()));
const LibUtilities::BasisKey B0(LibUtilities::eModified_A, order0,
LibUtilities::PointsKey(points0,LibUtilities::eGaussLobattoLegendre));
const LibUtilities::BasisKey B1(LibUtilities::eModified_B, order1,
LibUtilities::PointsKey(points1,LibUtilities::eGaussRadauMAlpha1Beta0));
m_xmap = MemoryManager<StdRegions::StdTriExp>::AllocateSharedPtr(B0,B1);
SetUpCoeffs(m_xmap->GetNcoeffs());
}
Nektar::SpatialDomains::TriGeom::TriGeom ( const int  id,
const SegGeomSharedPtr  edges[],
const StdRegions::Orientation  eorient[] 
)

Copy the edge shared pointers.

Definition at line 127 of file TriGeom.cpp.

References ASSERTL0, Nektar::StdRegions::eForwards, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eTriangle, Nektar::SpatialDomains::Geometry::GetBasis(), Nektar::SpatialDomains::Geometry::GetVertex(), kNedges, Nektar::SpatialDomains::Geometry::m_coordim, m_edges, m_eorient, m_fid, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_shapeType, m_verts, Nektar::SpatialDomains::Geometry::m_xmap, and Nektar::SpatialDomains::Geometry::SetUpCoeffs().

:
Geometry2D(edges[0]->GetVertex(0)->GetCoordim()),
m_fid(id)
{
/// Copy the edge shared pointers.
m_edges.insert(m_edges.begin(), edges, edges+TriGeom::kNedges);
for(int j=0; j <kNedges; ++j)
{
if(eorient[j] == StdRegions::eForwards)
{
m_verts.push_back(edges[j]->GetVertex(0));
}
else
{
m_verts.push_back(edges[j]->GetVertex(1));
}
}
for (int j=0; j<kNedges; ++j)
{
m_eorient[j] = eorient[j];
}
m_coordim = edges[0]->GetVertex(0)->GetCoordim();
ASSERTL0(m_coordim > 1,"Cannot call function with dim == 1");
int order0 = edges[0]->GetBasis(0)->GetNumModes();
int points0 = edges[0]->GetBasis(0)->GetNumPoints();
int order1 = max(order0,
max(edges[1]->GetBasis(0)->GetNumModes(),
edges[2]->GetBasis(0)->GetNumModes()));
int points1 = max(points0,
max(edges[1]->GetBasis(0)->GetNumPoints(),
edges[2]->GetBasis(0)->GetNumPoints()));
const LibUtilities::BasisKey B0(LibUtilities::eModified_A, order0,
LibUtilities::PointsKey(points0,LibUtilities::eGaussLobattoLegendre));
const LibUtilities::BasisKey B1(LibUtilities::eModified_B, order1,
LibUtilities::PointsKey(points1,LibUtilities::eGaussRadauMAlpha1Beta0));
m_xmap = MemoryManager<StdRegions::StdTriExp>::AllocateSharedPtr(B0,B1);
SetUpCoeffs(m_xmap->GetNcoeffs());
}
Nektar::SpatialDomains::TriGeom::TriGeom ( const int  id,
const SegGeomSharedPtr  edges[],
const StdRegions::Orientation  eorient[],
const CurveSharedPtr curve 
)

Copy the edge shared pointers.

Definition at line 180 of file TriGeom.cpp.

References ASSERTL0, Nektar::StdRegions::eForwards, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eModified_B, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eTriangle, Nektar::SpatialDomains::Geometry::GetBasis(), Nektar::LibUtilities::BasisKey::GetPointsKey(), Nektar::SpatialDomains::Geometry::GetVertex(), Nektar::SpatialDomains::Geometry::GetXmap(), Nektar::LibUtilities::Interp2D(), kNedges, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, m_edges, m_eorient, m_fid, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_shapeType, m_verts, Nektar::SpatialDomains::Geometry::m_xmap, npts, Nektar::LibUtilities::PointsManager(), and Nektar::SpatialDomains::Geometry::SetUpCoeffs().

:
Geometry2D(edges[0]->GetVertex(0)->GetCoordim()),
m_fid(id)
{
/// Copy the edge shared pointers.
m_edges.insert(m_edges.begin(), edges, edges+TriGeom::kNedges);
for(int j=0; j <kNedges; ++j)
{
if(eorient[j] == StdRegions::eForwards)
{
m_verts.push_back(edges[j]->GetVertex(0));
}
else
{
m_verts.push_back(edges[j]->GetVertex(1));
}
}
for (int j=0; j<kNedges; ++j)
{
m_eorient[j] = eorient[j];
}
m_coordim = edges[0]->GetVertex(0)->GetCoordim();
ASSERTL0(m_coordim > 1,"Cannot call function with dim == 1");
int order0 = edges[0]->GetBasis(0)->GetNumModes();
int points0 = edges[0]->GetBasis(0)->GetNumPoints();
int order1 = max(order0,
max(edges[1]->GetBasis(0)->GetNumModes(),
edges[2]->GetBasis(0)->GetNumModes()));
int points1 = max(points0,
max(edges[1]->GetBasis(0)->GetNumPoints(),
edges[2]->GetBasis(0)->GetNumPoints()));
const LibUtilities::BasisKey B0(LibUtilities::eModified_A, order0,
LibUtilities::PointsKey(points0,LibUtilities::eGaussLobattoLegendre));
const LibUtilities::BasisKey B1(LibUtilities::eModified_B, order1,
LibUtilities::PointsKey(points1,LibUtilities::eGaussRadauMAlpha1Beta0));
m_xmap = MemoryManager<StdRegions::StdTriExp>::AllocateSharedPtr(B0,B1);
SetUpCoeffs(m_xmap->GetNcoeffs());
for(int i = 0; i < m_coordim; ++i)
{
LibUtilities::PointsKey(2, curve->m_ptype)]->GetPointsDim();
// Deal with 2D points type separately (e.g. electrostatic or
// Fekete points).
if (pdim == 2)
{
int N = curve->m_points.size();
int nEdgePts = (-1+(int)sqrt(static_cast<NekDouble>(8*N+1)))/2;
ASSERTL0(nEdgePts*(nEdgePts+1)/2 == N,
"NUMPOINTS must be a triangle number for 2D basis.");
for (int j = 0; j < kNedges; ++j)
{
ASSERTL0(edges[j]->GetXmap()->GetNcoeffs() == nEdgePts,
"Number of edge points does not correspond "
"to number of face points.");
}
// Create a StdNodalTriExp.
const LibUtilities::PointsKey P0(
const LibUtilities::PointsKey P1(
const LibUtilities::BasisKey T0(
LibUtilities::eOrtho_A, nEdgePts, P0);
const LibUtilities::BasisKey T1(
LibUtilities::eOrtho_B, nEdgePts, P1);
MemoryManager<StdRegions::StdNodalTriExp>::AllocateSharedPtr(
T0, T1, curve->m_ptype);
Array<OneD, NekDouble> phys(max(t->GetTotPoints(),
m_xmap->GetTotPoints()));
for (int j = 0; j < N; ++j)
{
phys[j] = (curve->m_points[j]->GetPtr())[i];
}
Array<OneD, NekDouble> tmp(nEdgePts*nEdgePts);
t->BwdTrans(phys, tmp);
// Interpolate points to standard region.
B0.GetPointsKey(),B1.GetPointsKey(),
phys);
// Forwards transform to get coefficient space.
m_xmap->FwdTrans(phys, m_coeffs[i]);
}
else if (pdim == 1)
{
int npts = curve->m_points.size();
int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
Array<OneD,NekDouble> tmp(npts);
LibUtilities::PointsKey curveKey(nEdgePts, curve->m_ptype);
// Sanity checks:
// - Curved faces should have square number of points;
// - Each edge should have sqrt(npts) points.
ASSERTL0(nEdgePts*nEdgePts == npts,
"NUMPOINTS should be a square number");
for (int j = 0; j < kNedges; ++j)
{
ASSERTL0(edges[j]->GetXmap()->GetNcoeffs() == nEdgePts,
"Number of edge points does not correspond "
"to number of face points.");
}
for (int j = 0; j < npts; ++j)
{
tmp[j] = (curve->m_points[j]->GetPtr())[i];
}
// Interpolate curve points to standard triangle points.
Array<OneD, NekDouble> phys(m_xmap->GetTotPoints());
LibUtilities::Interp2D(curveKey,curveKey,tmp,
B0.GetPointsKey(),B1.GetPointsKey(),
phys);
// Forwards transform to get coefficient space.
m_xmap->FwdTrans(phys, m_coeffs[i]);
}
else
{
ASSERTL0(false, "Only 1D/2D points distributions supported.");
}
}
}
Nektar::SpatialDomains::TriGeom::TriGeom ( const TriGeom in)

Definition at line 329 of file TriGeom.cpp.

References kNedges, m_edges, m_elmtMap, m_eorient, m_fid, Nektar::SpatialDomains::Geometry::m_globalID, m_ownData, m_ownVerts, Nektar::SpatialDomains::Geometry::m_shapeType, and m_verts.

{
// From Geomtry
m_shapeType = in.m_shapeType;
// From TriFaceComponent
m_fid = in.m_fid;
m_ownVerts = in.m_ownVerts;
std::list<CompToElmt>::const_iterator def;
for(def = in.m_elmtMap.begin(); def != in.m_elmtMap.end(); def++)
{
m_elmtMap.push_back(*def);
}
// From TriGeom
m_verts = in.m_verts;
m_edges = in.m_edges;
for (int i = 0; i < kNedges; i++)
{
m_eorient[i] = in.m_eorient[i];
}
m_ownData = in.m_ownData;
}
Nektar::SpatialDomains::TriGeom::~TriGeom ( )

Definition at line 358 of file TriGeom.cpp.

{
}

Member Function Documentation

NekDouble Nektar::SpatialDomains::TriGeom::GetCoord ( const int  i,
const Array< OneD, const NekDouble > &  Lcoord 
)

Given local collapsed coordinate Lcoord return the value of physical coordinate in direction i

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 367 of file TriGeom.cpp.

References ASSERTL1, Nektar::SpatialDomains::ePtsFilled, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_state, and Nektar::SpatialDomains::Geometry::m_xmap.

Referenced by v_GetCoord().

{
"Geometry is not in physical space");
Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
m_xmap->BwdTrans(m_coeffs[i], tmp);
return m_xmap->PhysEvaluate(Lcoord, tmp);
}
StdRegions::Orientation Nektar::SpatialDomains::TriGeom::GetFaceOrientation ( const TriGeom face1,
const TriGeom face2 
)
static

Definition at line 379 of file TriGeom.cpp.

References m_verts.

Referenced by Nektar::MultiRegions::DisContField3D::FindPeriodicFaces().

{
return GetFaceOrientation(face1.m_verts, face2.m_verts);
}
StdRegions::Orientation Nektar::SpatialDomains::TriGeom::GetFaceOrientation ( const PointGeomVector face1,
const PointGeomVector face2 
)
static

Definition at line 386 of file TriGeom.cpp.

References ASSERTL0, Nektar::StdRegions::eDir1BwdDir1_Dir2BwdDir2, Nektar::StdRegions::eDir1BwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1BwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1FwdDir1_Dir2FwdDir2, Nektar::StdRegions::eDir1FwdDir2_Dir2BwdDir1, Nektar::StdRegions::eDir1FwdDir2_Dir2FwdDir1, and Nektar::StdRegions::eNoOrientation.

{
int i, j, vmap[3] = {-1,-1,-1};
NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
// For periodic faces, we calculate the vector between the centre
// points of the two faces. (For connected faces this will be
// zero). We can then use this to determine alignment later in the
// algorithm.
for (i = 0; i < 3; ++i)
{
cx += (*face2[i])(0) - (*face1[i])(0);
cy += (*face2[i])(1) - (*face1[i])(1);
cz += (*face2[i])(2) - (*face1[i])(2);
}
cx /= 3;
cy /= 3;
cz /= 3;
// Now construct a mapping which takes us from the vertices of one
// face to the other. That is, vertex j of face2 corresponds to
// vertex vmap[j] of face1.
for (i = 0; i < 3; ++i)
{
x = (*face1[i])(0);
y = (*face1[i])(1);
z = (*face1[i])(2);
for (j = 0; j < 3; ++j)
{
x1 = (*face2[j])(0)-cx;
y1 = (*face2[j])(1)-cy;
z1 = (*face2[j])(2)-cz;
if (sqrt((x1-x)*(x1-x)+(y1-y)*(y1-y)+(z1-z)*(z1-z)) < 1e-8)
{
vmap[j] = i;
break;
}
}
}
if (vmap[1] == (vmap[0]+1) % 3)
{
switch (vmap[0])
{
case 0:
break;
case 1:
break;
case 2:
break;
}
}
else
{
switch (vmap[0])
{
case 0:
break;
case 1:
break;
case 2:
break;
}
}
ASSERTL0(false, "Unable to determine triangle orientation");
}
void Nektar::SpatialDomains::TriGeom::v_AddElmtConnected ( int  gvo_id,
int  locid 
)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 467 of file TriGeom.cpp.

References m_elmtMap.

{
CompToElmt ee(gvo_id,locid);
m_elmtMap.push_back(ee);
}
bool Nektar::SpatialDomains::TriGeom::v_ContainsPoint ( const Array< OneD, const NekDouble > &  gloCoord,
NekDouble  tol = 0.0 
)
protectedvirtual

Determines if a point specified in global coordinates is located within this tetrahedral geometry.

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 817 of file TriGeom.cpp.

References Nektar::SpatialDomains::Geometry::GetCoordim().

Referenced by v_ContainsPoint().

{
Array<OneD,NekDouble> locCoord(GetCoordim(),0.0);
return v_ContainsPoint(gloCoord,locCoord,tol);
}
bool Nektar::SpatialDomains::TriGeom::v_ContainsPoint ( const Array< OneD, const NekDouble > &  gloCoord,
Array< OneD, NekDouble > &  locCoord,
NekDouble  tol 
)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 828 of file TriGeom.cpp.

References v_ContainsPoint().

{
NekDouble resid;
return v_ContainsPoint(gloCoord,stdCoord,tol,resid);
}
bool Nektar::SpatialDomains::TriGeom::v_ContainsPoint ( const Array< OneD, const NekDouble > &  gloCoord,
Array< OneD, NekDouble > &  locCoord,
NekDouble  tol,
NekDouble resid 
)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 836 of file TriGeom.cpp.

References ASSERTL1, and Nektar::SpatialDomains::Geometry::GetLocCoords().

{
ASSERTL1(gloCoord.num_elements() >= 2,
"Two dimensional geometry expects at least two coordinates.");
resid = GetLocCoords(gloCoord, stdCoord);
if (stdCoord[0] >= -(1+tol) && stdCoord[1] >= -(1+tol)
&& stdCoord[0] + stdCoord[1] <= tol)
{
return true;
}
return false;
}
void Nektar::SpatialDomains::TriGeom::v_FillGeom ( )
protectedvirtual

Put all quadrature information into edge structure.

Note verts and edges are listed according to anticlockwise convention but points in _coeffs have to be in array format from left to right.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 594 of file TriGeom.cpp.

References Nektar::SpatialDomains::ePtsFilled, kNedges, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, m_edges, m_eorient, Nektar::SpatialDomains::Geometry::m_state, and Nektar::SpatialDomains::Geometry::m_xmap.

Referenced by v_GenGeomFactors(), and v_GetLocCoords().

{
// check to see if geometry structure is already filled
{
int i,j,k;
int nEdgeCoeffs = m_xmap->GetEdgeNcoeffs(0);
Array<OneD, unsigned int> mapArray (nEdgeCoeffs);
Array<OneD, int> signArray(nEdgeCoeffs);
for(i = 0; i < kNedges; i++)
{
m_edges[i]->FillGeom();
m_xmap->GetEdgeToElementMap(i,m_eorient[i],mapArray,signArray);
nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
for(j = 0 ; j < m_coordim; j++)
{
for(k = 0; k < nEdgeCoeffs; k++)
{
m_coeffs[j][mapArray[k]] =
signArray[k] * m_edges[i]->GetCoeffs(j)[k];
}
}
}
}
}
void Nektar::SpatialDomains::TriGeom::v_GenGeomFactors ( )
protectedvirtual

Set up GeoFac for this geometry using Coord quadrature distribution

Implements Nektar::SpatialDomains::Geometry.

Definition at line 555 of file TriGeom.cpp.

References Nektar::SpatialDomains::eDeformed, Nektar::SpatialDomains::ePtsFilled, Nektar::SpatialDomains::eRegular, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, Nektar::SpatialDomains::Geometry::m_geomFactors, Nektar::SpatialDomains::Geometry::m_geomFactorsState, Nektar::SpatialDomains::Geometry::m_xmap, and v_FillGeom().

{
{
GeomType Gtype = eRegular;
// check to see if expansions are linear
for(int i = 0; i < m_coordim; ++i)
{
if(m_xmap->GetBasisNumModes(0) != 2 ||
m_xmap->GetBasisNumModes(1) != 2)
{
Gtype = eDeformed;
}
}
m_geomFactors = MemoryManager<GeomFactors>::AllocateSharedPtr(
Gtype, m_coordim, m_xmap, m_coeffs);
}
}
const LibUtilities::BasisSharedPtr Nektar::SpatialDomains::TriGeom::v_GetBasis ( const int  i)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 519 of file TriGeom.cpp.

References Nektar::SpatialDomains::Geometry::m_xmap.

{
return m_xmap->GetBasis(i);
}
StdRegions::Orientation Nektar::SpatialDomains::TriGeom::v_GetCartesianEorient ( const int  i) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 751 of file TriGeom.cpp.

References ASSERTL2, Nektar::StdRegions::eBackwards, Nektar::StdRegions::eForwards, and m_eorient.

{
ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
if(i < 2)
{
return m_eorient[i];
}
else
{
{
}
else
{
}
}
}
NekDouble Nektar::SpatialDomains::TriGeom::v_GetCoord ( const int  i,
const Array< OneD, const NekDouble > &  Lcoord 
)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 546 of file TriGeom.cpp.

References GetCoord().

{
return GetCoord(i,Lcoord);
}
int Nektar::SpatialDomains::TriGeom::v_GetCoordim ( void  ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 510 of file TriGeom.cpp.

References Nektar::SpatialDomains::Geometry::m_coordim.

{
return m_coordim;
}
const Geometry1DSharedPtr Nektar::SpatialDomains::TriGeom::v_GetEdge ( int  i) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 731 of file TriGeom.cpp.

References ASSERTL2, and m_edges.

{
ASSERTL2((i >=0) && (i <= 2),"Edge id must be between 0 and 3");
return m_edges[i];
}
const LibUtilities::BasisSharedPtr Nektar::SpatialDomains::TriGeom::v_GetEdgeBasis ( const int  i)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 528 of file TriGeom.cpp.

References ASSERTL1, and Nektar::SpatialDomains::Geometry::m_xmap.

{
ASSERTL1(i <= 2, "edge is out of range");
if(i == 0)
{
return m_xmap->GetBasis(0);
}
else
{
return m_xmap->GetBasis(1);
}
}
int Nektar::SpatialDomains::TriGeom::v_GetEid ( int  i) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 701 of file TriGeom.cpp.

References ASSERTL2, and m_edges.

{
ASSERTL2((i >=0) && (i <= 2),"Edge id must be between 0 and 2");
return m_edges[i]->GetEid();
}
StdRegions::Orientation Nektar::SpatialDomains::TriGeom::v_GetEorient ( const int  i) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 741 of file TriGeom.cpp.

References ASSERTL2, and m_eorient.

{
ASSERTL2((i >=0) && (i <= 2),"Edge id must be between 0 and 2");
return m_eorient[i];
}
int Nektar::SpatialDomains::TriGeom::v_GetFid ( ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 501 of file TriGeom.cpp.

References m_fid.

{
return m_fid;
}
NekDouble Nektar::SpatialDomains::TriGeom::v_GetLocCoords ( const Array< OneD, const NekDouble > &  coords,
Array< OneD, NekDouble > &  Lcoords 
)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 630 of file TriGeom.cpp.

References Nektar::SpatialDomains::PointGeom::dot(), Nektar::SpatialDomains::eRegular, Nektar::SpatialDomains::Geometry::GetMetricInfo(), Vmath::Imin(), Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, m_verts, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::PointGeom::Mult(), Nektar::SpatialDomains::Geometry2D::NewtonIterationForLocCoord(), npts, Vmath::Sadd(), Nektar::SpatialDomains::PointGeom::Sub(), v_FillGeom(), Vmath::Vmul(), and Vmath::Vvtvp().

{
NekDouble resid = 0.0;
// calculate local coordinate for coord
if(GetMetricInfo()->GetGtype() == eRegular)
{
NekDouble coords2 = (m_coordim == 3)? coords[2]: 0.0;
PointGeom dv1, dv2, norm, orth1, orth2;
PointGeom xin(m_coordim,0,coords[0],coords[1],coords2);
// Calculate edge vectors from 0-1 and 0-2 edges.
dv1.Sub(*m_verts[1],*m_verts[0]);
dv2.Sub(*m_verts[2],*m_verts[0]);
// Obtain normal to plane in which dv1 and dv2 lie
norm.Mult(dv1,dv2);
// Obtain vector which are proportional to normal of dv1 and dv2.
orth1.Mult(norm,dv1);
orth2.Mult(norm,dv2);
// Start with vector of desired points minus vertex_0
xin -= *m_verts[0];
// Calculate length using L/|dv1| = (x-v0).n1/(dv1.n1) for coordiante 1
// Then rescale to [-1,1].
Lcoords[0] = xin.dot(orth2)/dv1.dot(orth2);
Lcoords[0] = 2*Lcoords[0]-1;
Lcoords[1] = xin.dot(orth1)/dv2.dot(orth1);
Lcoords[1] = 2*Lcoords[1]-1;
}
else
{
// Determine nearest point of coords to values in m_xmap
int npts = m_xmap->GetTotPoints();
Array<OneD, NekDouble> ptsx(npts), ptsy(npts);
Array<OneD, NekDouble> tmpx(npts), tmpy(npts);
m_xmap->BwdTrans(m_coeffs[0], ptsx);
m_xmap->BwdTrans(m_coeffs[1], ptsy);
const Array<OneD, const NekDouble> za = m_xmap->GetPoints(0);
const Array<OneD, const NekDouble> zb = m_xmap->GetPoints(1);
//guess the first local coords based on nearest point
Vmath::Sadd(npts, -coords[0], ptsx,1,tmpx,1);
Vmath::Sadd(npts, -coords[1], ptsy,1,tmpy,1);
Vmath::Vmul (npts, tmpx,1,tmpx,1,tmpx,1);
Vmath::Vvtvp(npts, tmpy,1,tmpy,1,tmpx,1,tmpx,1);
int min_i = Vmath::Imin(npts,tmpx,1);
Lcoords[0] = za[min_i%za.num_elements()];
Lcoords[1] = zb[min_i/za.num_elements()];
// recover cartesian coordinate from collapsed coordinate.
Lcoords[0] = (1.0+Lcoords[0])*(1.0-Lcoords[1])/2 -1.0;
// Perform newton iteration to find local coordinates
NewtonIterationForLocCoord(coords, ptsx, ptsy, Lcoords,resid);
}
return resid;
}
int Nektar::SpatialDomains::TriGeom::v_GetNumEdges ( ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 807 of file TriGeom.cpp.

References kNedges.

{
return kNedges;
}
int Nektar::SpatialDomains::TriGeom::v_GetNumVerts ( ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 798 of file TriGeom.cpp.

References kNverts.

{
return kNverts;
}
PointGeomSharedPtr Nektar::SpatialDomains::TriGeom::v_GetVertex ( int  i) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 721 of file TriGeom.cpp.

References ASSERTL2, and m_verts.

{
ASSERTL2((i >=0) && (i <= 2),"Vertex id must be between 0 and 2");
return m_verts[i];
}
int Nektar::SpatialDomains::TriGeom::v_GetVid ( int  i) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 711 of file TriGeom.cpp.

References ASSERTL2, and m_verts.

{
ASSERTL2((i >=0) && (i <= 2),"Vertex id must be between 0 and 2");
return m_verts[i]->GetVid();
}
bool Nektar::SpatialDomains::TriGeom::v_IsElmtConnected ( int  gvo_id,
int  locid 
) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 486 of file TriGeom.cpp.

References Nektar::StdRegions::find(), and m_elmtMap.

{
std::list<CompToElmt>::const_iterator def;
CompToElmt ee(gvo_id,locid);
def = find(m_elmtMap.begin(),m_elmtMap.end(),ee);
// Found the element connectivity object in the list
return (def != m_elmtMap.end());
}
int Nektar::SpatialDomains::TriGeom::v_NumElmtConnected ( ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 477 of file TriGeom.cpp.

References m_elmtMap.

{
return int(m_elmtMap.size());
}
void Nektar::SpatialDomains::TriGeom::v_SetOwnData ( )
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 583 of file TriGeom.cpp.

References m_ownData.

{
m_ownData = true;
}
int Nektar::SpatialDomains::TriGeom::v_WhichEdge ( SegGeomSharedPtr  edge)
protectedvirtual

Return the edge number of the given edge.

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 775 of file TriGeom.cpp.

References Nektar::iterator, and m_edges.

{
int returnval = -1;
int i;
for (i=0,edgeIter = m_edges.begin(); edgeIter != m_edges.end(); ++edgeIter,++i)
{
if (*edgeIter == edge)
{
returnval = i;
break;
}
}
return returnval;
}

Member Data Documentation

const int Nektar::SpatialDomains::TriGeom::kNedges = 3
static

Get the orientation of face1.

Definition at line 98 of file TriGeom.h.

Referenced by Nektar::SpatialDomains::MeshGraph2D::ReadElements(), Nektar::SpatialDomains::MeshGraph3D::ReadFaces(), TriGeom(), v_FillGeom(), and v_GetNumEdges().

const int Nektar::SpatialDomains::TriGeom::kNverts = 3
static

Definition at line 99 of file TriGeom.h.

Referenced by TriGeom(), and v_GetNumVerts().

SegGeomVector Nektar::SpatialDomains::TriGeom::m_edges
protected

Definition at line 112 of file TriGeom.h.

Referenced by TriGeom(), v_FillGeom(), v_GetEdge(), v_GetEid(), and v_WhichEdge().

std::list<CompToElmt> Nektar::SpatialDomains::TriGeom::m_elmtMap
protected

Definition at line 116 of file TriGeom.h.

Referenced by TriGeom(), v_AddElmtConnected(), v_IsElmtConnected(), and v_NumElmtConnected().

StdRegions::Orientation Nektar::SpatialDomains::TriGeom::m_eorient[kNedges]
protected

Definition at line 113 of file TriGeom.h.

Referenced by TriGeom(), v_FillGeom(), v_GetCartesianEorient(), and v_GetEorient().

int Nektar::SpatialDomains::TriGeom::m_fid
protected

Definition at line 114 of file TriGeom.h.

Referenced by TriGeom(), and v_GetFid().

bool Nektar::SpatialDomains::TriGeom::m_ownData
private

Definition at line 193 of file TriGeom.h.

Referenced by TriGeom(), and v_SetOwnData().

bool Nektar::SpatialDomains::TriGeom::m_ownVerts
protected

Definition at line 115 of file TriGeom.h.

Referenced by TriGeom().

PointGeomVector Nektar::SpatialDomains::TriGeom::m_verts
protected

Definition at line 111 of file TriGeom.h.

Referenced by GetFaceOrientation(), TriGeom(), v_GetLocCoords(), v_GetVertex(), and v_GetVid().