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 Member Functions | 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 (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)
 
- 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 More...
 
int GetEdgeFaceMap (int i, int j) const
 
void FillGeom ()
 Put all quadrature information into face/edge structure and backward transform. More...
 
NekDouble GetLocCoords (const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
 
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. More...
 
void SetOwnData ()
 
const LibUtilities::BasisSharedPtr GetBasis (const int i)
 Return the j-th basis of the i-th co-ordinate dimension. More...
 
const LibUtilities::PointsKeyVector GetPointsKeys ()
 
void Reset (CurveMap &curvedEdges, CurveMap &curvedFaces)
 

Static Public Member Functions

static StdRegions::Orientation GetFaceOrientation (const QuadGeom &face1, const QuadGeom &face2)
 Get the orientation of face1. More...
 
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. More...
 
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. More...
 
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)
 
virtual void v_Reset (CurveMap &curvedEdges, CurveMap &curvedFaces)
 Reset this geometry object: unset the current state and remove allocated GeomFactors. More...
 
- 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
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)
 Initialise the m_coeffs array. More...
 

Protected Attributes

PointGeomVector m_verts
 
SegGeomVector m_edges
 
StdRegions::Orientation m_eorient [kNedges]
 
int m_fid
 
bool m_ownVerts
 
std::list< CompToElmtm_elmtMap
 
CurveSharedPtr m_curve
 
- Protected Attributes inherited from Nektar::SpatialDomains::Geometry
int m_coordim
 coordinate dimension More...
 
GeomFactorsSharedPtr m_geomFactors
 
GeomState m_geomFactorsState
 
StdRegions::StdExpansionSharedPtr m_xmap
 
GeomState m_state
 enum identifier to determine if quad points are filled More...
 
GeomType m_geomType
 
LibUtilities::ShapeType m_shapeType
 
int m_globalID
 
Array< OneD, Array< OneD,
NekDouble > > 
m_coeffs
 

Private Member Functions

void SetUpXmap ()
 

Private Attributes

bool m_ownData
 Boolean indicating whether object owns the data. More...
 

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 ( 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 62 of file QuadGeom.cpp.

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

65  :
66  Geometry2D(verts[0]->GetCoordim()),
67  m_fid(id)
68  {
69  m_globalID = id;
71 
72  /// Copy the vert shared pointers.
73  m_verts.insert(m_verts.begin(), verts, verts+QuadGeom::kNverts);
74 
75  /// Copy the edge shared pointers.
76  m_edges.insert(m_edges.begin(), edges, edges+QuadGeom::kNedges);
77 
78 
79  for (int j=0; j<kNedges; ++j)
80  {
81  m_eorient[j] = eorient[j];
82  }
83 
84  m_coordim = verts[0]->GetCoordim();
85  ASSERTL0(m_coordim > 1,
86  "Cannot call function with dim == 1");
87 
88  SetUpXmap();
89  SetUpCoeffs(m_xmap->GetNcoeffs());
90  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
LibUtilities::ShapeType m_shapeType
Definition: Geometry.h:177
void SetUpCoeffs(const int nCoeffs)
Initialise the m_coeffs array.
Definition: Geometry.h:484
StdRegions::Orientation m_eorient[kNedges]
Definition: QuadGeom.h:107
int m_coordim
coordinate dimension
Definition: Geometry.h:169
Nektar::SpatialDomains::QuadGeom::QuadGeom ( const int  id,
const SegGeomSharedPtr  edges[],
const StdRegions::Orientation  eorient[] 
)

Copy the edge shared pointers.

Definition at line 142 of file QuadGeom.cpp.

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

144  :
145  Geometry2D(edges[0]->GetVertex(0)->GetCoordim()),
146  m_fid(id)
147  {
148  int j;
149 
150  m_globalID = m_fid;
152 
153  /// Copy the edge shared pointers.
154  m_edges.insert(m_edges.begin(), edges, edges+QuadGeom::kNedges);
155 
156  for(j=0; j <kNedges; ++j)
157  {
158  if(eorient[j] == StdRegions::eForwards)
159  {
160  m_verts.push_back(edges[j]->GetVertex(0));
161  }
162  else
163  {
164  m_verts.push_back(edges[j]->GetVertex(1));
165  }
166  }
167 
168  for (j=0; j<kNedges; ++j)
169  {
170  m_eorient[j] = eorient[j];
171  }
172 
173  m_coordim = edges[0]->GetVertex(0)->GetCoordim();
174  ASSERTL0(m_coordim > 1,
175  "Cannot call function with dim == 1");
176 
177  SetUpXmap();
178  SetUpCoeffs(m_xmap->GetNcoeffs());
179  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
PointGeomSharedPtr GetVertex(int i) const
Definition: Geometry.h:348
LibUtilities::ShapeType m_shapeType
Definition: Geometry.h:177
void SetUpCoeffs(const int nCoeffs)
Initialise the m_coeffs array.
Definition: Geometry.h:484
StdRegions::Orientation m_eorient[kNedges]
Definition: QuadGeom.h:107
int m_coordim
coordinate dimension
Definition: Geometry.h:169
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 96 of file QuadGeom.cpp.

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

99  :
100  Geometry2D(edges[0]->GetVertex(0)->GetCoordim()),
101  m_fid(id),
102  m_curve(curve)
103  {
104  int j;
105 
106  m_globalID = m_fid;
107 
109 
110  /// Copy the edge shared pointers.
111  m_edges.insert(m_edges.begin(), edges, edges+QuadGeom::kNedges);
112 
113  for(j=0; j <kNedges; ++j)
114  {
115  if(eorient[j] == StdRegions::eForwards)
116  {
117  m_verts.push_back(edges[j]->GetVertex(0));
118  }
119  else
120  {
121  m_verts.push_back(edges[j]->GetVertex(1));
122  }
123  }
124 
125  for (j=0; j<kNedges; ++j)
126  {
127  m_eorient[j] = eorient[j];
128  }
129 
130  m_coordim = edges[0]->GetVertex(0)->GetCoordim();
131  ASSERTL0(m_coordim > 1,
132  "Cannot call function with dim == 1");
133 
134  SetUpXmap();
135  SetUpCoeffs(m_xmap->GetNcoeffs());
136  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
PointGeomSharedPtr GetVertex(int i) const
Definition: Geometry.h:348
LibUtilities::ShapeType m_shapeType
Definition: Geometry.h:177
void SetUpCoeffs(const int nCoeffs)
Initialise the m_coeffs array.
Definition: Geometry.h:484
StdRegions::Orientation m_eorient[kNedges]
Definition: QuadGeom.h:107
int m_coordim
coordinate dimension
Definition: Geometry.h:169
Nektar::SpatialDomains::QuadGeom::QuadGeom ( const QuadGeom in)

Definition at line 185 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.

186  {
187  // From Geometry
188  m_shapeType = in.m_shapeType;
189 
190  // From QuadFaceComponent
191  m_fid = in.m_fid;
192  m_globalID = m_fid;
193  m_ownVerts = in.m_ownVerts;
194  std::list<CompToElmt>::const_iterator def;
195  for(def = in.m_elmtMap.begin(); def != in.m_elmtMap.end(); def++)
196  {
197  m_elmtMap.push_back(*def);
198  }
199 
200  // From QuadGeom
201  m_verts = in.m_verts;
202  m_edges = in.m_edges;
203  for (int i = 0; i < kNedges; i++)
204  {
205  m_eorient[i] = in.m_eorient[i];
206  }
207  m_ownData = in.m_ownData;
208  }
std::list< CompToElmt > m_elmtMap
Definition: QuadGeom.h:110
LibUtilities::ShapeType m_shapeType
Definition: Geometry.h:177
StdRegions::Orientation m_eorient[kNedges]
Definition: QuadGeom.h:107
bool m_ownData
Boolean indicating whether object owns the data.
Definition: QuadGeom.h:193
Nektar::SpatialDomains::QuadGeom::~QuadGeom ( )

Definition at line 214 of file QuadGeom.cpp.

215  {
216  }

Member Function Documentation

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

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

243  {
245  "Geometry is not in physical space");
246 
247  Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
248  m_xmap->BwdTrans(m_coeffs[i], tmp);
249 
250  return m_xmap->PhysEvaluate(Lcoord, tmp);
251  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Definition: Geometry.h:180
Geometric information has been generated.
GeomState m_state
enum identifier to determine if quad points are filled
Definition: Geometry.h:175
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
StdRegions::Orientation Nektar::SpatialDomains::QuadGeom::GetFaceOrientation ( const QuadGeom face1,
const QuadGeom face2 
)
static

Get the orientation of face1.

Definition at line 253 of file QuadGeom.cpp.

References m_verts.

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

256  {
257  return GetFaceOrientation(face1.m_verts, face2.m_verts);
258  }
static StdRegions::Orientation GetFaceOrientation(const QuadGeom &face1, const QuadGeom &face2)
Get the orientation of face1.
Definition: QuadGeom.cpp:253
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 264 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.

267  {
268  int i, j, vmap[4] = {-1,-1,-1,-1};
269  NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
270 
271  // For periodic faces, we calculate the vector between the centre
272  // points of the two faces. (For connected faces this will be
273  // zero). We can then use this to determine alignment later in the
274  // algorithm.
275  for (i = 0; i < 4; ++i)
276  {
277  cx += (*face2[i])(0) - (*face1[i])(0);
278  cy += (*face2[i])(1) - (*face1[i])(1);
279  cz += (*face2[i])(2) - (*face1[i])(2);
280  }
281  cx /= 4;
282  cy /= 4;
283  cz /= 4;
284 
285  // Now construct a mapping which takes us from the vertices of one
286  // face to the other. That is, vertex j of face2 corresponds to
287  // vertex vmap[j] of face1.
288  for (i = 0; i < 4; ++i)
289  {
290  x = (*face1[i])(0);
291  y = (*face1[i])(1);
292  z = (*face1[i])(2);
293  for (j = 0; j < 4; ++j)
294  {
295  x1 = (*face2[j])(0)-cx;
296  y1 = (*face2[j])(1)-cy;
297  z1 = (*face2[j])(2)-cz;
298  if (sqrt((x1-x)*(x1-x)+(y1-y)*(y1-y)+(z1-z)*(z1-z)) < 1e-8)
299  {
300  vmap[j] = i;
301  break;
302  }
303  }
304  }
305 
306  // Use the mapping to determine the eight alignment options between
307  // faces.
308  if (vmap[1] == (vmap[0]+1) % 4)
309  {
310  switch (vmap[0])
311  {
312  case 0:
314  break;
315  case 1:
317  break;
318  case 2:
320  break;
321  case 3:
323  break;
324  }
325  }
326  else
327  {
328  switch (vmap[0])
329  {
330  case 0:
332  break;
333  case 1:
335  break;
336  case 2:
338  break;
339  case 3:
341  break;
342  }
343  }
344 
345  ASSERTL0(false, "unable to determine face orientation");
347  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
double NekDouble
void Nektar::SpatialDomains::QuadGeom::SetUpXmap ( )
private

Definition at line 218 of file QuadGeom.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eModified_A, Nektar::SpatialDomains::Geometry::GetBasis(), m_edges, and Nektar::SpatialDomains::Geometry::m_xmap.

Referenced by QuadGeom(), and v_Reset().

219  {
220  int order0 = max(m_edges[0]->GetBasis(0)->GetNumModes(),
221  m_edges[2]->GetBasis(0)->GetNumModes());
222  int order1 = max(m_edges[1]->GetBasis(0)->GetNumModes(),
223  m_edges[3]->GetBasis(0)->GetNumModes());
224 
225  const LibUtilities::BasisKey B0(
227  LibUtilities::PointsKey(
229  const LibUtilities::BasisKey B1(
231  LibUtilities::PointsKey(
233 
235  ::AllocateSharedPtr(B0,B1);
236  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
Principle Modified Functions .
Definition: BasisType.h:49
const LibUtilities::BasisSharedPtr GetBasis(const int i)
Return the j-th basis of the i-th co-ordinate dimension.
Definition: Geometry.h:475
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50
void Nektar::SpatialDomains::QuadGeom::v_AddElmtConnected ( int  gvo_id,
int  locid 
)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 353 of file QuadGeom.cpp.

References m_elmtMap.

354  {
355  CompToElmt ee(gvo_id,locid);
356  m_elmtMap.push_back(ee);
357  }
std::list< CompToElmt > m_elmtMap
Definition: QuadGeom.h:110
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 788 of file QuadGeom.cpp.

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

Referenced by v_ContainsPoint().

790  {
791  Array<OneD,NekDouble> locCoord(GetCoordim(),0.0);
792  return v_ContainsPoint(gloCoord,locCoord,tol);
793 
794  }
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
Definition: QuadGeom.cpp:788
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 798 of file QuadGeom.cpp.

References v_ContainsPoint().

801  {
802  NekDouble resid;
803  return v_ContainsPoint(gloCoord,stdCoord,tol,resid);
804  }
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
Definition: QuadGeom.cpp:788
double NekDouble
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 806 of file QuadGeom.cpp.

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

810  {
811  ASSERTL1(gloCoord.num_elements() >= 2,
812  "Two dimensional geometry expects at least two coordinates.");
813 
814  resid = GetLocCoords(gloCoord, stdCoord);
815  if (stdCoord[0] >= -(1+tol) && stdCoord[1] >= -(1+tol)
816  && stdCoord[0] <= (1+tol) && stdCoord[1] <= (1+tol))
817  {
818  return true;
819  }
820  return false;
821  }
NekDouble GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Definition: Geometry.h:450
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
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 527 of file QuadGeom.cpp.

References ASSERTL0, Nektar::SpatialDomains::ePtsFilled, Nektar::SpatialDomains::Geometry::GetXmap(), Nektar::LibUtilities::Interp2D(), kNedges, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, m_curve, m_edges, m_eorient, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_state, Nektar::SpatialDomains::Geometry::m_xmap, and npts.

Referenced by v_GenGeomFactors(), and v_GetLocCoords().

528  {
529  // check to see if geometry structure is already filled
530  if(m_state != ePtsFilled)
531  {
532  int i,j,k;
533  int nEdgeCoeffs;
534 
535  if (m_curve)
536  {
537  int npts = m_curve->m_points.size();
538  int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
539  Array<OneD, NekDouble> tmp(npts);
540  Array<OneD, NekDouble> tmp2(m_xmap->GetTotPoints());
541  LibUtilities::PointsKey curveKey(
542  nEdgePts, m_curve->m_ptype);
543 
544  // Sanity checks:
545  // - Curved faces should have square number of points;
546  // - Each edge should have sqrt(npts) points.
547  ASSERTL0(nEdgePts * nEdgePts == npts,
548  "NUMPOINTS should be a square number in"
549  " quadrilteral "
550  + boost::lexical_cast<string>(m_globalID));
551 
552  for (i = 0; i < kNedges; ++i)
553  {
554  ASSERTL0(
555  m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
556  "Number of edge points does not correspond to "
557  "number of face points in quadrilateral "
558  + boost::lexical_cast<string>(m_globalID));
559  }
560 
561  for (i = 0; i < m_coordim; ++i)
562  {
563  for (j = 0; j < npts; ++j)
564  {
565  tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
566  }
567 
568  // Interpolate m_curve points to GLL points
570  curveKey, curveKey, tmp,
571  m_xmap->GetBasis(0)->GetPointsKey(),
572  m_xmap->GetBasis(1)->GetPointsKey(),
573  tmp2);
574 
575  // Forwards transform to get coefficient space.
576  m_xmap->FwdTrans(tmp2, m_coeffs[i]);
577  }
578  }
579 
580  // Now fill in edges.
581  Array<OneD, unsigned int> mapArray;
582  Array<OneD, int> signArray;
583 
584  for(i = 0; i < kNedges; i++)
585  {
586  m_edges[i]->FillGeom();
587  m_xmap->GetEdgeToElementMap(i,m_eorient[i],
588  mapArray,signArray);
589 
590  nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
591 
592  for(j = 0; j < m_coordim; j++)
593  {
594  for(k = 0; k < nEdgeCoeffs; k++)
595  {
596  m_coeffs[j][mapArray[k]]
597  = signArray[k]*(m_edges[i]->GetCoeffs(j))[k];
598  }
599  }
600  }
601 
603  }
604  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
StdRegions::StdExpansionSharedPtr GetXmap() const
Definition: Geometry.h:383
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis...
Definition: Interp.cpp:116
static std::string npts
Definition: InputFld.cpp:43
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Definition: Geometry.h:180
Geometric information has been generated.
StdRegions::Orientation m_eorient[kNedges]
Definition: QuadGeom.h:107
GeomState m_state
enum identifier to determine if quad points are filled
Definition: Geometry.h:175
int m_coordim
coordinate dimension
Definition: Geometry.h:169
void Nektar::SpatialDomains::QuadGeom::v_GenGeomFactors ( )
protectedvirtual

Set up GeoFac for this geometry using Coord quadrature distribution

Implements Nektar::SpatialDomains::Geometry.

Definition at line 445 of file QuadGeom.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), 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().

446  {
448  {
449  int i;
450  GeomType Gtype = eRegular;
451 
453 
454  // We will first check whether we have a regular or deformed
455  // geometry. We will define regular as those cases where the
456  // Jacobian and the metric terms of the derivative are constants
457  // (i.e. not coordinate dependent)
458 
459  // Check to see if expansions are linear
460  // If not linear => deformed geometry
461  for(i = 0; i < m_coordim; ++i)
462  {
463  if((m_xmap->GetBasisNumModes(0) != 2)||
464  (m_xmap->GetBasisNumModes(1) != 2))
465  {
466  Gtype = eDeformed;
467  }
468  }
469 
470  // For linear expansions, the mapping from standard to local
471  // element is given by the relation:
472  // x_i = 0.25 * [ ( x_i^A + x_i^B + x_i^C + x_i^D) +
473  // (-x_i^A + x_i^B + x_i^C - x_i^D)*xi_1 +
474  // (-x_i^A - x_i^B + x_i^C + x_i^D)*xi_2 +
475  // ( x_i^A - x_i^B + x_i^C - x_i^D)*xi_1*xi_2 ]
476  //
477  // The jacobian of the transformation and the metric terms
478  // dxi_i/dx_j, involve only terms of the form dx_i/dxi_j (both
479  // for coordim == 2 or 3). Inspecting the formula above, it can
480  // be appreciated that the derivatives dx_i/dxi_j will be
481  // constant, if the coefficient of the non-linear term is zero.
482  //
483  // That is why for regular geometry, we require
484  //
485  // x_i^A - x_i^B + x_i^C - x_i^D = 0
486  //
487  // or equivalently
488  //
489  // x_i^A - x_i^B = x_i^D - x_i^C
490  //
491  // This corresponds to quadrilaterals which are paralellograms.
492  if(Gtype == eRegular)
493  {
494  for(i = 0; i < m_coordim; i++)
495  {
496  if( fabs( (*m_verts[0])(i) - (*m_verts[1])(i) +
497  (*m_verts[2])(i) - (*m_verts[3])(i) )
499  {
500  Gtype = eDeformed;
501  break;
502  }
503  }
504  }
505 
507  Gtype, m_coordim, m_xmap, m_coeffs);
509  }
510  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
GeomFactorsSharedPtr m_geomFactors
Definition: Geometry.h:170
static const NekDouble kNekZeroTol
virtual void v_FillGeom()
Put all quadrature information into edge structure.
Definition: QuadGeom.cpp:527
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Definition: Geometry.h:180
Geometry is straight-sided with constant geometric factors.
Geometric information has been generated.
GeomType
Indicates the type of element geometry.
Geometry is curved or has non-constant factors.
int m_coordim
coordinate dimension
Definition: Geometry.h:169
const LibUtilities::BasisSharedPtr Nektar::SpatialDomains::QuadGeom::v_GetBasis ( const int  i)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 410 of file QuadGeom.cpp.

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

411  {
412  return m_xmap->GetBasis(i);
413  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
StdRegions::Orientation Nektar::SpatialDomains::QuadGeom::v_GetCartesianEorient ( const int  i) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 726 of file QuadGeom.cpp.

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

727  {
728  ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
729  if(i < 2)
730  {
731  return m_eorient[i];
732  }
733  else
734  {
736  {
737  return StdRegions::eBackwards;
738  }
739  else
740  {
741  return StdRegions::eForwards;
742  }
743  }
744  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:213
StdRegions::Orientation m_eorient[kNedges]
Definition: QuadGeom.h:107
NekDouble Nektar::SpatialDomains::QuadGeom::v_GetCoord ( const int  i,
const Array< OneD, const NekDouble > &  Lcoord 
)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 436 of file QuadGeom.cpp.

References GetCoord().

437  {
438  return GetCoord(i,Lcoord);
439  }
NekDouble GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord)
Definition: QuadGeom.cpp:241
int Nektar::SpatialDomains::QuadGeom::v_GetCoordim ( void  ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 401 of file QuadGeom.cpp.

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

402  {
403  return m_coordim;
404  }
int m_coordim
coordinate dimension
Definition: Geometry.h:169
const Geometry1DSharedPtr Nektar::SpatialDomains::QuadGeom::v_GetEdge ( int  i) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 706 of file QuadGeom.cpp.

References ASSERTL2, and m_edges.

707  {
708  ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
709  return m_edges[i];
710  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:213
const LibUtilities::BasisSharedPtr Nektar::SpatialDomains::QuadGeom::v_GetEdgeBasis ( const int  i)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 419 of file QuadGeom.cpp.

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

420  {
421  ASSERTL1(i <= 3,"edge is out of range");
422 
423  if (i == 0 || i == 2)
424  {
425  return m_xmap->GetBasis(0);
426  }
427  else
428  {
429  return m_xmap->GetBasis(1);
430  }
431  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
int Nektar::SpatialDomains::QuadGeom::v_GetEid ( int  i) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 676 of file QuadGeom.cpp.

References ASSERTL2, and m_edges.

677  {
678  ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
679  return m_edges[i]->GetEid();
680  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:213
StdRegions::Orientation Nektar::SpatialDomains::QuadGeom::v_GetEorient ( const int  i) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 716 of file QuadGeom.cpp.

References ASSERTL2, and m_eorient.

717  {
718  ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
719  return m_eorient[i];
720  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:213
StdRegions::Orientation m_eorient[kNedges]
Definition: QuadGeom.h:107
int Nektar::SpatialDomains::QuadGeom::v_GetFid ( ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 392 of file QuadGeom.cpp.

References m_fid.

393  {
394  return m_fid;
395  }
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 609 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().

611  {
612  NekDouble resid = 0.0;
613  if(GetMetricInfo()->GetGtype() == eRegular)
614  {
615  NekDouble coords2 = (m_coordim == 3)? coords[2]: 0.0;
616  PointGeom dv1, dv2, norm, orth1, orth2;
617  PointGeom xin(m_coordim,0,coords[0],coords[1],coords2);
618 
619  // Calculate edge vectors from 0-1 and 0-3 edges.
620  dv1.Sub(*m_verts[1],*m_verts[0]);
621  dv2.Sub(*m_verts[3],*m_verts[0]);
622 
623  // Obtain normal to plane in which dv1 and dv2 lie
624  norm.Mult(dv1,dv2);
625 
626  // Obtain vector which are normal to dv1 and dv2.
627  orth1.Mult(norm,dv1);
628  orth2.Mult(norm,dv2);
629 
630  // Start with vector of desired points minus vertex_0
631  xin -= *m_verts[0];
632 
633  // Calculate length using L/|dv1| = (x-v0).n1/(dv1.n1) for coordiante 1
634  // Then rescale to [-1,1].
635  Lcoords[0] = xin.dot(orth2)/dv1.dot(orth2);
636  Lcoords[0] = 2*Lcoords[0]-1;
637  Lcoords[1] = xin.dot(orth1)/dv2.dot(orth1);
638  Lcoords[1] = 2*Lcoords[1]-1;
639  }
640  else
641  {
643 
644  // Determine nearest point of coords to values in m_xmap
645  int npts = m_xmap->GetTotPoints();
646  Array<OneD, NekDouble> ptsx(npts), ptsy(npts);
647  Array<OneD, NekDouble> tmpx(npts), tmpy(npts);
648 
649  m_xmap->BwdTrans(m_coeffs[0], ptsx);
650  m_xmap->BwdTrans(m_coeffs[1], ptsy);
651 
652  const Array<OneD, const NekDouble> za = m_xmap->GetPoints(0);
653  const Array<OneD, const NekDouble> zb = m_xmap->GetPoints(1);
654 
655  //guess the first local coords based on nearest point
656  Vmath::Sadd(npts, -coords[0], ptsx,1,tmpx,1);
657  Vmath::Sadd(npts, -coords[1], ptsy,1,tmpy,1);
658  Vmath::Vmul (npts, tmpx,1,tmpx,1,tmpx,1);
659  Vmath::Vvtvp(npts, tmpy,1,tmpy,1,tmpx,1,tmpx,1);
660 
661  int min_i = Vmath::Imin(npts,tmpx,1);
662 
663  Lcoords[0] = za[min_i%za.num_elements()];
664  Lcoords[1] = zb[min_i/za.num_elements()];
665 
666  // Perform newton iteration to find local coordinates
667  NewtonIterationForLocCoord(coords, ptsx, ptsy, Lcoords,resid);
668  }
669  return resid;
670  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
int Imin(int n, const T *x, const int incx)
Return the index of the minimum element in x.
Definition: Vmath.cpp:833
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:428
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)
Definition: Geometry2D.cpp:100
static std::string npts
Definition: InputFld.cpp:43
double NekDouble
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
Definition: Vmath.cpp:301
virtual void v_FillGeom()
Put all quadrature information into edge structure.
Definition: QuadGeom.cpp:527
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Definition: Geometry.h:180
Geometry is straight-sided with constant geometric factors.
GeomFactorsSharedPtr GetMetricInfo()
Definition: Geometry.h:299
int m_coordim
coordinate dimension
Definition: Geometry.h:169
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:169
int Nektar::SpatialDomains::QuadGeom::v_GetNumEdges ( ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 782 of file QuadGeom.cpp.

References kNedges.

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 773 of file QuadGeom.cpp.

References kNverts.

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

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 696 of file QuadGeom.cpp.

References ASSERTL2, and m_verts.

697  {
698  ASSERTL2((i >=0) && (i <= 3),"Vertex id must be between 0 and 3");
699  return m_verts[i];
700  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:213
int Nektar::SpatialDomains::QuadGeom::v_GetVid ( int  i) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 686 of file QuadGeom.cpp.

References ASSERTL2, and m_verts.

687  {
688  ASSERTL2((i >=0) && (i <= 3),"Verted id must be between 0 and 3");
689  return m_verts[i]->GetVid();
690  }
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:213
bool Nektar::SpatialDomains::QuadGeom::v_IsElmtConnected ( int  gvo_id,
int  locid 
) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 372 of file QuadGeom.cpp.

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

373  {
374  std::list<CompToElmt>::const_iterator def;
375  CompToElmt ee(gvo_id,locid);
376 
377  def = find(m_elmtMap.begin(),m_elmtMap.end(),ee);
378 
379  // Found the element connectivity object in the list
380  if(def != m_elmtMap.end())
381  {
382  return(true);
383  }
384 
385  return(false);
386  }
std::list< CompToElmt > m_elmtMap
Definition: QuadGeom.h:110
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:315
int Nektar::SpatialDomains::QuadGeom::v_NumElmtConnected ( ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 363 of file QuadGeom.cpp.

References m_elmtMap.

364  {
365  return int(m_elmtMap.size());
366  }
std::list< CompToElmt > m_elmtMap
Definition: QuadGeom.h:110
void Nektar::SpatialDomains::QuadGeom::v_Reset ( CurveMap curvedEdges,
CurveMap curvedFaces 
)
protectedvirtual

Reset this geometry object: unset the current state and remove allocated GeomFactors.

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 823 of file QuadGeom.cpp.

References Nektar::iterator, m_curve, m_edges, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::Geometry::SetUpCoeffs(), SetUpXmap(), and Nektar::SpatialDomains::Geometry::v_Reset().

825  {
826  Geometry::v_Reset(curvedEdges, curvedFaces);
827  CurveMap::iterator it = curvedFaces.find(m_globalID);
828 
829  if (it != curvedFaces.end())
830  {
831  m_curve = it->second;
832  }
833 
834  for (int i = 0; i < 4; ++i)
835  {
836  m_edges[i]->Reset(curvedEdges, curvedFaces);
837  }
838 
839  SetUpXmap();
840  SetUpCoeffs(m_xmap->GetNcoeffs());
841  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state and remove allocated GeomFactors.
Definition: Geometry.cpp:307
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
void SetUpCoeffs(const int nCoeffs)
Initialise the m_coeffs array.
Definition: Geometry.h:484
void Nektar::SpatialDomains::QuadGeom::v_SetOwnData ( )
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 516 of file QuadGeom.cpp.

References m_ownData.

517  {
518  m_ownData = true;
519  }
bool m_ownData
Boolean indicating whether object owns the data.
Definition: QuadGeom.h:193
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 750 of file QuadGeom.cpp.

References Nektar::iterator, and m_edges.

751  {
752  int returnval = -1;
753 
754  SegGeomVector::iterator edgeIter;
755  int i;
756 
757  for (i=0,edgeIter = m_edges.begin(); edgeIter != m_edges.end(); ++edgeIter,++i)
758  {
759  if (*edgeIter == edge)
760  {
761  returnval = i;
762  break;
763  }
764  }
765 
766  return returnval;
767  }
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator

Member Data Documentation

const int Nektar::SpatialDomains::QuadGeom::kNedges = 4
static
const int Nektar::SpatialDomains::QuadGeom::kNverts = 4
static
CurveSharedPtr Nektar::SpatialDomains::QuadGeom::m_curve
protected

Definition at line 111 of file QuadGeom.h.

Referenced by v_FillGeom(), and v_Reset().

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

Definition at line 106 of file QuadGeom.h.

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

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

Definition at line 110 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 107 of file QuadGeom.h.

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

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

Definition at line 108 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 193 of file QuadGeom.h.

Referenced by QuadGeom(), and v_SetOwnData().

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

Definition at line 109 of file QuadGeom.h.

Referenced by QuadGeom().

PointGeomVector Nektar::SpatialDomains::QuadGeom::m_verts
protected
const std::string Nektar::SpatialDomains::QuadGeom::XMLElementType
static

Definition at line 102 of file QuadGeom.h.