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::QuadGeom Class Reference

#include <QuadGeom.h>

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

Public Member Functions

 QuadGeom ()
 QuadGeom (int id, const int coordim)
 QuadGeom (const int id, const PointGeomSharedPtr verts[], const SegGeomSharedPtr edges[], const StdRegions::Orientation eorient[])
 QuadGeom (const int id, const SegGeomSharedPtr edges[], const StdRegions::Orientation eorient[])
 QuadGeom (const int id, const SegGeomSharedPtr edges[], const StdRegions::Orientation eorient[], const CurveSharedPtr &curve)
 QuadGeom (const QuadGeom &in)
 ~QuadGeom ()
NekDouble GetCoord (const int i, const Array< OneD, const NekDouble > &Lcoord)
 Given local collapsed coordinate Lcoord return the value of physical coordinate in direction i.
- 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 GetNumVerts () const
PointGeomSharedPtr GetVertex (int i) const
StdRegions::Orientation GetEorient (const int i) const
StdRegions::Orientation GetPorient (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 QuadGeom &face1, const QuadGeom &face2)
 Get the orientation of face1.
static StdRegions::Orientation GetFaceOrientation (const PointGeomVector &face1, const PointGeomVector &face2)

Static Public Attributes

static const int kNverts = 4
static const int kNedges = 4
static const std::string XMLElementType

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)
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)
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 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
 Boolean indicating whether object owns the data.

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 60 of file QuadGeom.h.

Constructor & Destructor Documentation

Nektar::SpatialDomains::QuadGeom::QuadGeom ( )
Nektar::SpatialDomains::QuadGeom::QuadGeom ( int  id,
const int  coordim 
)

Definition at line 62 of file QuadGeom.cpp.

References Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eModified_A, m_fid, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_xmap, and Nektar::SpatialDomains::Geometry::SetUpCoeffs().

:
Geometry2D(coordim)
{
const LibUtilities::BasisKey B(LibUtilities::eModified_A, 2,
LibUtilities::PointsKey(3,LibUtilities::eGaussLobattoLegendre));
m_globalID = m_fid = id;
m_xmap = MemoryManager<StdRegions::StdQuadExp>::AllocateSharedPtr(B,B);
SetUpCoeffs(m_xmap->GetNcoeffs());
}
Nektar::SpatialDomains::QuadGeom::QuadGeom ( 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 78 of file QuadGeom.cpp.

References ASSERTL0, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eQuadrilateral, Nektar::SpatialDomains::Geometry::GetBasis(), kNedges, kNverts, Nektar::SpatialDomains::Geometry::m_coordim, m_edges, m_eorient, 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)
{
m_globalID = id;
/// Copy the vert shared pointers.
m_verts.insert(m_verts.begin(), verts, verts+QuadGeom::kNverts);
/// Copy the edge shared pointers.
m_edges.insert(m_edges.begin(), edges, edges+QuadGeom::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 = max(edges[0]->GetBasis(0)->GetNumModes(),
edges[2]->GetBasis(0)->GetNumModes());
int points0 = max(edges[0]->GetBasis(0)->GetNumPoints(),
edges[2]->GetBasis(0)->GetNumPoints());
int order1 = max(edges[1]->GetBasis(0)->GetNumModes(),
edges[3]->GetBasis(0)->GetNumModes());
int points1 = max(edges[1]->GetBasis(0)->GetNumPoints(),
edges[3]->GetBasis(0)->GetNumPoints());
const LibUtilities::BasisKey B0(LibUtilities::eModified_A, order0,
LibUtilities::PointsKey(points0,LibUtilities::eGaussLobattoLegendre));
const LibUtilities::BasisKey B1(LibUtilities::eModified_A, order1,
LibUtilities::PointsKey(points1,LibUtilities::eGaussLobattoLegendre));
m_xmap = MemoryManager<StdRegions::StdQuadExp>::AllocateSharedPtr(B0,B1);
SetUpCoeffs(m_xmap->GetNcoeffs());
}
Nektar::SpatialDomains::QuadGeom::QuadGeom ( const int  id,
const SegGeomSharedPtr  edges[],
const StdRegions::Orientation  eorient[] 
)

Copy the edge shared pointers.

Definition at line 220 of file QuadGeom.cpp.

References ASSERTL0, Nektar::StdRegions::eForwards, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eQuadrilateral, 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)
{
int j;
/// Copy the edge shared pointers.
m_edges.insert(m_edges.begin(), edges, edges+QuadGeom::kNedges);
for(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 (j=0; j<kNedges; ++j)
{
m_eorient[j] = eorient[j];
}
m_coordim = edges[0]->GetVertex(0)->GetCoordim();
"Cannot call function with dim == 1");
int order0 = max(edges[0]->GetBasis(0)->GetNumModes(),
edges[2]->GetBasis(0)->GetNumModes());
int points0 = max(edges[0]->GetBasis(0)->GetNumPoints(),
edges[2]->GetBasis(0)->GetNumPoints());
int order1 = max(edges[1]->GetBasis(0)->GetNumModes(),
edges[3]->GetBasis(0)->GetNumModes());
int points1 = max(edges[1]->GetBasis(0)->GetNumPoints(),
edges[3]->GetBasis(0)->GetNumPoints());
const LibUtilities::BasisKey B0(LibUtilities::eModified_A, order0,
LibUtilities::PointsKey(points0,LibUtilities::eGaussLobattoLegendre));
const LibUtilities::BasisKey B1(LibUtilities::eModified_A, order1,
LibUtilities::PointsKey(points1,LibUtilities::eGaussLobattoLegendre));
m_xmap = MemoryManager<StdRegions::StdQuadExp>::AllocateSharedPtr(B0,B1);
SetUpCoeffs(m_xmap->GetNcoeffs());
}
Nektar::SpatialDomains::QuadGeom::QuadGeom ( const int  id,
const SegGeomSharedPtr  edges[],
const StdRegions::Orientation  eorient[],
const CurveSharedPtr curve 
)

Copy the edge shared pointers.

Definition at line 126 of file QuadGeom.cpp.

References ASSERTL0, Nektar::StdRegions::eForwards, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eQuadrilateral, 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, and Nektar::SpatialDomains::Geometry::SetUpCoeffs().

:
Geometry2D(edges[0]->GetVertex(0)->GetCoordim()),
m_fid(id)
{
int j;
/// Copy the edge shared pointers.
m_edges.insert(m_edges.begin(), edges, edges+QuadGeom::kNedges);
for(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 (j=0; j<kNedges; ++j)
{
m_eorient[j] = eorient[j];
}
m_coordim = edges[0]->GetVertex(0)->GetCoordim();
"Cannot call function with dim == 1");
int order0 = max(edges[0]->GetBasis(0)->GetNumModes(),
edges[2]->GetBasis(0)->GetNumModes());
int points0 = max(edges[0]->GetBasis(0)->GetNumPoints(),
edges[2]->GetBasis(0)->GetNumPoints());
int order1 = max(edges[1]->GetBasis(0)->GetNumModes(),
edges[3]->GetBasis(0)->GetNumModes());
int points1 = max(edges[1]->GetBasis(0)->GetNumPoints(),
edges[3]->GetBasis(0)->GetNumPoints());
const LibUtilities::BasisKey B0(LibUtilities::eModified_A, order0,
LibUtilities::PointsKey(points0,LibUtilities::eGaussLobattoLegendre));
const LibUtilities::BasisKey B1(LibUtilities::eModified_A, order1,
LibUtilities::PointsKey(points1,LibUtilities::eGaussLobattoLegendre));
m_xmap = MemoryManager<StdRegions::StdQuadExp>::AllocateSharedPtr(B0,B1);
SetUpCoeffs(m_xmap->GetNcoeffs());
for(int i = 0; i < m_coordim; ++i)
{
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 (j = 0; j < kNedges; ++j)
{
ASSERTL0(edges[j]->GetXmap()->GetNcoeffs() == nEdgePts,
"Number of edge points does not correspond "
"to number of face points.");
}
for (j = 0; j < npts; ++j)
{
tmp[j] = (curve->m_points[j]->GetPtr())[i];
}
// Interpolate curve points to GLL points
Array<OneD, NekDouble> tmp2(points0*points1);
LibUtilities::Interp2D(curveKey,curveKey,tmp,
B0.GetPointsKey(),B1.GetPointsKey(),
tmp2);
// Forwards transform to get coefficient space.
m_xmap->FwdTrans(tmp2, m_coeffs[i]);
}
}
Nektar::SpatialDomains::QuadGeom::QuadGeom ( const QuadGeom in)

Definition at line 277 of file QuadGeom.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 Geometry
m_shapeType = in.m_shapeType;
// From QuadFaceComponent
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 QuadGeom
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::QuadGeom::~QuadGeom ( )

Definition at line 306 of file QuadGeom.cpp.

{
}

Member Function Documentation

NekDouble Nektar::SpatialDomains::QuadGeom::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 314 of file QuadGeom.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::QuadGeom::GetFaceOrientation ( const QuadGeom face1,
const QuadGeom face2 
)
static

Get the orientation of face1.

Definition at line 326 of file QuadGeom.cpp.

References m_verts.

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

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

Calculate the orientation of face2 to face1 (note this is not face1 to face2!).

Definition at line 337 of file QuadGeom.cpp.

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

{
int i, j, vmap[4] = {-1,-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 < 4; ++i)
{
cx += (*face2[i])(0) - (*face1[i])(0);
cy += (*face2[i])(1) - (*face1[i])(1);
cz += (*face2[i])(2) - (*face1[i])(2);
}
cx /= 4;
cy /= 4;
cz /= 4;
// 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 < 4; ++i)
{
x = (*face1[i])(0);
y = (*face1[i])(1);
z = (*face1[i])(2);
for (j = 0; j < 4; ++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;
}
}
}
// Use the mapping to determine the eight alignment options between
// faces.
if (vmap[1] == (vmap[0]+1) % 4)
{
switch (vmap[0])
{
case 0:
break;
case 1:
break;
case 2:
break;
case 3:
break;
}
}
else
{
switch (vmap[0])
{
case 0:
break;
case 1:
break;
case 2:
break;
case 3:
break;
}
}
ASSERTL0(false, "unable to determine face orientation");
}
void Nektar::SpatialDomains::QuadGeom::v_AddElmtConnected ( int  gvo_id,
int  locid 
)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 426 of file QuadGeom.cpp.

References m_elmtMap.

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

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 815 of file QuadGeom.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::QuadGeom::v_ContainsPoint ( const Array< OneD, const NekDouble > &  gloCoord,
Array< OneD, NekDouble > &  locCoord,
NekDouble  tol 
)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 825 of file QuadGeom.cpp.

References v_ContainsPoint().

{
NekDouble resid;
return v_ContainsPoint(gloCoord,stdCoord,tol,resid);
}
bool Nektar::SpatialDomains::QuadGeom::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 833 of file QuadGeom.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] <= (1+tol) && stdCoord[1] <= (1+tol))
{
return true;
}
return false;
}
void Nektar::SpatialDomains::QuadGeom::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 600 of file QuadGeom.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;
Array<OneD, unsigned int> mapArray;
Array<OneD, int> signArray;
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::QuadGeom::v_GenGeomFactors ( )
protectedvirtual

Set up GeoFac for this geometry using Coord quadrature distribution

Implements Nektar::SpatialDomains::Geometry.

Definition at line 518 of file QuadGeom.cpp.

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

{
{
int i;
GeomType Gtype = eRegular;
// We will first check whether we have a regular or deformed
// geometry. We will define regular as those cases where the
// Jacobian and the metric terms of the derivative are constants
// (i.e. not coordinate dependent)
// Check to see if expansions are linear
// If not linear => deformed geometry
for(i = 0; i < m_coordim; ++i)
{
if((m_xmap->GetBasisNumModes(0) != 2)||
(m_xmap->GetBasisNumModes(1) != 2))
{
Gtype = eDeformed;
}
}
// For linear expansions, the mapping from standard to local
// element is given by the relation:
// x_i = 0.25 * [ ( x_i^A + x_i^B + x_i^C + x_i^D) +
// (-x_i^A + x_i^B + x_i^C - x_i^D)*xi_1 +
// (-x_i^A - x_i^B + x_i^C + x_i^D)*xi_2 +
// ( x_i^A - x_i^B + x_i^C - x_i^D)*xi_1*xi_2 ]
//
// The jacobian of the transformation and the metric terms
// dxi_i/dx_j, involve only terms of the form dx_i/dxi_j (both
// for coordim == 2 or 3). Inspecting the formula above, it can
// be appreciated that the derivatives dx_i/dxi_j will be
// constant, if the coefficient of the non-linear term is zero.
//
// That is why for regular geometry, we require
//
// x_i^A - x_i^B + x_i^C - x_i^D = 0
//
// or equivalently
//
// x_i^A - x_i^B = x_i^D - x_i^C
//
// This corresponds to quadrilaterals which are paralellograms.
if(Gtype == eRegular)
{
for(i = 0; i < m_coordim; i++)
{
if( fabs( (*m_verts[0])(i) - (*m_verts[1])(i) +
(*m_verts[2])(i) - (*m_verts[3])(i) )
{
Gtype = eDeformed;
break;
}
}
}
m_geomFactors = MemoryManager<GeomFactors>::AllocateSharedPtr(
Gtype, m_coordim, m_xmap, m_coeffs);
}
}
const LibUtilities::BasisSharedPtr Nektar::SpatialDomains::QuadGeom::v_GetBasis ( const int  i)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 483 of file QuadGeom.cpp.

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

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

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 753 of file QuadGeom.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::QuadGeom::v_GetCoord ( const int  i,
const Array< OneD, const NekDouble > &  Lcoord 
)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 509 of file QuadGeom.cpp.

References GetCoord().

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 474 of file QuadGeom.cpp.

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

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

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 733 of file QuadGeom.cpp.

References ASSERTL2, and m_edges.

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

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 492 of file QuadGeom.cpp.

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

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

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 703 of file QuadGeom.cpp.

References ASSERTL2, and m_edges.

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

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 743 of file QuadGeom.cpp.

References ASSERTL2, and m_eorient.

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

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 465 of file QuadGeom.cpp.

References m_fid.

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 636 of file QuadGeom.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;
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-3 edges.
dv1.Sub(*m_verts[1],*m_verts[0]);
dv2.Sub(*m_verts[3],*m_verts[0]);
// Obtain normal to plane in which dv1 and dv2 lie
norm.Mult(dv1,dv2);
// Obtain vector which are normal to 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()];
// Perform newton iteration to find local coordinates
NewtonIterationForLocCoord(coords, ptsx, ptsy, Lcoords,resid);
}
return resid;
}
int Nektar::SpatialDomains::QuadGeom::v_GetNumEdges ( ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 809 of file QuadGeom.cpp.

References kNedges.

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 800 of file QuadGeom.cpp.

References kNverts.

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

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 723 of file QuadGeom.cpp.

References ASSERTL2, and m_verts.

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 713 of file QuadGeom.cpp.

References ASSERTL2, and m_verts.

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 445 of file QuadGeom.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
if(def != m_elmtMap.end())
{
return(true);
}
return(false);
}
int Nektar::SpatialDomains::QuadGeom::v_NumElmtConnected ( ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 436 of file QuadGeom.cpp.

References m_elmtMap.

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 589 of file QuadGeom.cpp.

References m_ownData.

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

Return the edge number of the given edge.

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 777 of file QuadGeom.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::QuadGeom::kNedges = 4
static

Definition at line 103 of file QuadGeom.h.

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

const int Nektar::SpatialDomains::QuadGeom::kNverts = 4
static

Definition at line 102 of file QuadGeom.h.

Referenced by QuadGeom(), and v_GetNumVerts().

SegGeomVector Nektar::SpatialDomains::QuadGeom::m_edges
protected

Definition at line 108 of file QuadGeom.h.

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

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

Definition at line 112 of file QuadGeom.h.

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

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

Definition at line 109 of file QuadGeom.h.

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

int Nektar::SpatialDomains::QuadGeom::m_fid
protected

Definition at line 110 of file QuadGeom.h.

Referenced by QuadGeom(), and v_GetFid().

bool Nektar::SpatialDomains::QuadGeom::m_ownData
private

Boolean indicating whether object owns the data.

Definition at line 189 of file QuadGeom.h.

Referenced by QuadGeom(), and v_SetOwnData().

bool Nektar::SpatialDomains::QuadGeom::m_ownVerts
protected

Definition at line 111 of file QuadGeom.h.

Referenced by QuadGeom().

PointGeomVector Nektar::SpatialDomains::QuadGeom::m_verts
protected

Definition at line 107 of file QuadGeom.h.

Referenced by GetFaceOrientation(), QuadGeom(), v_GenGeomFactors(), v_GetLocCoords(), v_GetVertex(), and v_GetVid().

const std::string Nektar::SpatialDomains::QuadGeom::XMLElementType
static

Definition at line 104 of file QuadGeom.h.