Nektar++
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 245 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().

247  {
249  "Geometry is not in physical space");
250 
251  Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
252  m_xmap->BwdTrans(m_coeffs[i], tmp);
253 
254  return m_xmap->PhysEvaluate(Lcoord, tmp);
255  }
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 257 of file QuadGeom.cpp.

References m_verts.

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

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

271  {
272  int i, j, vmap[4] = {-1,-1,-1,-1};
273  NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
274 
275  // For periodic faces, we calculate the vector between the centre
276  // points of the two faces. (For connected faces this will be
277  // zero). We can then use this to determine alignment later in the
278  // algorithm.
279  for (i = 0; i < 4; ++i)
280  {
281  cx += (*face2[i])(0) - (*face1[i])(0);
282  cy += (*face2[i])(1) - (*face1[i])(1);
283  cz += (*face2[i])(2) - (*face1[i])(2);
284  }
285  cx /= 4;
286  cy /= 4;
287  cz /= 4;
288 
289  // Now construct a mapping which takes us from the vertices of one
290  // face to the other. That is, vertex j of face2 corresponds to
291  // vertex vmap[j] of face1.
292  for (i = 0; i < 4; ++i)
293  {
294  x = (*face1[i])(0);
295  y = (*face1[i])(1);
296  z = (*face1[i])(2);
297  for (j = 0; j < 4; ++j)
298  {
299  x1 = (*face2[j])(0)-cx;
300  y1 = (*face2[j])(1)-cy;
301  z1 = (*face2[j])(2)-cz;
302  if (sqrt((x1-x)*(x1-x)+(y1-y)*(y1-y)+(z1-z)*(z1-z)) < 1e-8)
303  {
304  vmap[j] = i;
305  break;
306  }
307  }
308  }
309 
310  // Use the mapping to determine the eight alignment options between
311  // faces.
312  if (vmap[1] == (vmap[0]+1) % 4)
313  {
314  switch (vmap[0])
315  {
316  case 0:
318  break;
319  case 1:
321  break;
322  case 2:
324  break;
325  case 3:
327  break;
328  }
329  }
330  else
331  {
332  switch (vmap[0])
333  {
334  case 0:
336  break;
337  case 1:
339  break;
340  case 2:
342  break;
343  case 3:
345  break;
346  }
347  }
348 
349  ASSERTL0(false, "unable to determine face orientation");
351  }
#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 points0 = max(m_edges[0]->GetBasis(0)->GetNumPoints(),
223  m_edges[2]->GetBasis(0)->GetNumPoints());
224  int order1 = max(m_edges[1]->GetBasis(0)->GetNumModes(),
225  m_edges[3]->GetBasis(0)->GetNumModes());
226  int points1 = max(m_edges[1]->GetBasis(0)->GetNumPoints(),
227  m_edges[3]->GetBasis(0)->GetNumPoints());
228 
229  const LibUtilities::BasisKey B0(
231  LibUtilities::PointsKey(
233  const LibUtilities::BasisKey B1(
235  LibUtilities::PointsKey(
237 
239  ::AllocateSharedPtr(B0,B1);
240  }
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 357 of file QuadGeom.cpp.

References m_elmtMap.

358  {
359  CompToElmt ee(gvo_id,locid);
360  m_elmtMap.push_back(ee);
361  }
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 792 of file QuadGeom.cpp.

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

Referenced by v_ContainsPoint().

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

References v_ContainsPoint().

805  {
806  NekDouble resid;
807  return v_ContainsPoint(gloCoord,stdCoord,tol,resid);
808  }
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
Definition: QuadGeom.cpp:792
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 810 of file QuadGeom.cpp.

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

814  {
815  ASSERTL1(gloCoord.num_elements() >= 2,
816  "Two dimensional geometry expects at least two coordinates.");
817 
818  resid = GetLocCoords(gloCoord, stdCoord);
819  if (stdCoord[0] >= -(1+tol) && stdCoord[1] >= -(1+tol)
820  && stdCoord[0] <= (1+tol) && stdCoord[1] <= (1+tol))
821  {
822  return true;
823  }
824  return false;
825  }
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 531 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().

532  {
533  // check to see if geometry structure is already filled
534  if(m_state != ePtsFilled)
535  {
536  int i,j,k;
537  int nEdgeCoeffs;
538 
539  if (m_curve)
540  {
541  int npts = m_curve->m_points.size();
542  int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
543  Array<OneD, NekDouble> tmp(npts);
544  Array<OneD, NekDouble> tmp2(m_xmap->GetTotPoints());
545  LibUtilities::PointsKey curveKey(
546  nEdgePts, m_curve->m_ptype);
547 
548  // Sanity checks:
549  // - Curved faces should have square number of points;
550  // - Each edge should have sqrt(npts) points.
551  ASSERTL0(nEdgePts * nEdgePts == npts,
552  "NUMPOINTS should be a square number in"
553  " quadrilteral "
554  + boost::lexical_cast<string>(m_globalID));
555 
556  for (i = 0; i < kNedges; ++i)
557  {
558  ASSERTL0(
559  m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
560  "Number of edge points does not correspond to "
561  "number of face points in quadrilateral "
562  + boost::lexical_cast<string>(m_globalID));
563  }
564 
565  for (i = 0; i < m_coordim; ++i)
566  {
567  for (j = 0; j < npts; ++j)
568  {
569  tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
570  }
571 
572  // Interpolate m_curve points to GLL points
574  curveKey, curveKey, tmp,
575  m_xmap->GetBasis(0)->GetPointsKey(),
576  m_xmap->GetBasis(1)->GetPointsKey(),
577  tmp2);
578 
579  // Forwards transform to get coefficient space.
580  m_xmap->FwdTrans(tmp2, m_coeffs[i]);
581  }
582  }
583 
584  // Now fill in edges.
585  Array<OneD, unsigned int> mapArray;
586  Array<OneD, int> signArray;
587 
588  for(i = 0; i < kNedges; i++)
589  {
590  m_edges[i]->FillGeom();
591  m_xmap->GetEdgeToElementMap(i,m_eorient[i],
592  mapArray,signArray);
593 
594  nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
595 
596  for(j = 0; j < m_coordim; j++)
597  {
598  for(k = 0; k < nEdgeCoeffs; k++)
599  {
600  m_coeffs[j][mapArray[k]]
601  = signArray[k]*(m_edges[i]->GetCoeffs(j))[k];
602  }
603  }
604  }
605 
607  }
608  }
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 449 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().

450  {
452  {
453  int i;
454  GeomType Gtype = eRegular;
455 
457 
458  // We will first check whether we have a regular or deformed
459  // geometry. We will define regular as those cases where the
460  // Jacobian and the metric terms of the derivative are constants
461  // (i.e. not coordinate dependent)
462 
463  // Check to see if expansions are linear
464  // If not linear => deformed geometry
465  for(i = 0; i < m_coordim; ++i)
466  {
467  if((m_xmap->GetBasisNumModes(0) != 2)||
468  (m_xmap->GetBasisNumModes(1) != 2))
469  {
470  Gtype = eDeformed;
471  }
472  }
473 
474  // For linear expansions, the mapping from standard to local
475  // element is given by the relation:
476  // x_i = 0.25 * [ ( x_i^A + x_i^B + x_i^C + x_i^D) +
477  // (-x_i^A + x_i^B + x_i^C - x_i^D)*xi_1 +
478  // (-x_i^A - x_i^B + x_i^C + x_i^D)*xi_2 +
479  // ( x_i^A - x_i^B + x_i^C - x_i^D)*xi_1*xi_2 ]
480  //
481  // The jacobian of the transformation and the metric terms
482  // dxi_i/dx_j, involve only terms of the form dx_i/dxi_j (both
483  // for coordim == 2 or 3). Inspecting the formula above, it can
484  // be appreciated that the derivatives dx_i/dxi_j will be
485  // constant, if the coefficient of the non-linear term is zero.
486  //
487  // That is why for regular geometry, we require
488  //
489  // x_i^A - x_i^B + x_i^C - x_i^D = 0
490  //
491  // or equivalently
492  //
493  // x_i^A - x_i^B = x_i^D - x_i^C
494  //
495  // This corresponds to quadrilaterals which are paralellograms.
496  if(Gtype == eRegular)
497  {
498  for(i = 0; i < m_coordim; i++)
499  {
500  if( fabs( (*m_verts[0])(i) - (*m_verts[1])(i) +
501  (*m_verts[2])(i) - (*m_verts[3])(i) )
503  {
504  Gtype = eDeformed;
505  break;
506  }
507  }
508  }
509 
511  Gtype, m_coordim, m_xmap, m_coeffs);
513  }
514  }
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:531
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 414 of file QuadGeom.cpp.

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

415  {
416  return m_xmap->GetBasis(i);
417  }
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 730 of file QuadGeom.cpp.

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

731  {
732  ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
733  if(i < 2)
734  {
735  return m_eorient[i];
736  }
737  else
738  {
740  {
741  return StdRegions::eBackwards;
742  }
743  else
744  {
745  return StdRegions::eForwards;
746  }
747  }
748  }
#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 440 of file QuadGeom.cpp.

References GetCoord().

441  {
442  return GetCoord(i,Lcoord);
443  }
NekDouble GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord)
Definition: QuadGeom.cpp:245
int Nektar::SpatialDomains::QuadGeom::v_GetCoordim ( void  ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 405 of file QuadGeom.cpp.

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

406  {
407  return m_coordim;
408  }
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 710 of file QuadGeom.cpp.

References ASSERTL2, and m_edges.

711  {
712  ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
713  return m_edges[i];
714  }
#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 423 of file QuadGeom.cpp.

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

424  {
425  ASSERTL1(i <= 3,"edge is out of range");
426 
427  if (i == 0 || i == 2)
428  {
429  return m_xmap->GetBasis(0);
430  }
431  else
432  {
433  return m_xmap->GetBasis(1);
434  }
435  }
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 680 of file QuadGeom.cpp.

References ASSERTL2, and m_edges.

681  {
682  ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
683  return m_edges[i]->GetEid();
684  }
#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 720 of file QuadGeom.cpp.

References ASSERTL2, and m_eorient.

721  {
722  ASSERTL2((i >=0) && (i <= 3),"Edge id must be between 0 and 3");
723  return m_eorient[i];
724  }
#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 396 of file QuadGeom.cpp.

References m_fid.

397  {
398  return m_fid;
399  }
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 613 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().

615  {
616  NekDouble resid = 0.0;
617  if(GetMetricInfo()->GetGtype() == eRegular)
618  {
619  NekDouble coords2 = (m_coordim == 3)? coords[2]: 0.0;
620  PointGeom dv1, dv2, norm, orth1, orth2;
621  PointGeom xin(m_coordim,0,coords[0],coords[1],coords2);
622 
623  // Calculate edge vectors from 0-1 and 0-3 edges.
624  dv1.Sub(*m_verts[1],*m_verts[0]);
625  dv2.Sub(*m_verts[3],*m_verts[0]);
626 
627  // Obtain normal to plane in which dv1 and dv2 lie
628  norm.Mult(dv1,dv2);
629 
630  // Obtain vector which are normal to dv1 and dv2.
631  orth1.Mult(norm,dv1);
632  orth2.Mult(norm,dv2);
633 
634  // Start with vector of desired points minus vertex_0
635  xin -= *m_verts[0];
636 
637  // Calculate length using L/|dv1| = (x-v0).n1/(dv1.n1) for coordiante 1
638  // Then rescale to [-1,1].
639  Lcoords[0] = xin.dot(orth2)/dv1.dot(orth2);
640  Lcoords[0] = 2*Lcoords[0]-1;
641  Lcoords[1] = xin.dot(orth1)/dv2.dot(orth1);
642  Lcoords[1] = 2*Lcoords[1]-1;
643  }
644  else
645  {
647 
648  // Determine nearest point of coords to values in m_xmap
649  int npts = m_xmap->GetTotPoints();
650  Array<OneD, NekDouble> ptsx(npts), ptsy(npts);
651  Array<OneD, NekDouble> tmpx(npts), tmpy(npts);
652 
653  m_xmap->BwdTrans(m_coeffs[0], ptsx);
654  m_xmap->BwdTrans(m_coeffs[1], ptsy);
655 
656  const Array<OneD, const NekDouble> za = m_xmap->GetPoints(0);
657  const Array<OneD, const NekDouble> zb = m_xmap->GetPoints(1);
658 
659  //guess the first local coords based on nearest point
660  Vmath::Sadd(npts, -coords[0], ptsx,1,tmpx,1);
661  Vmath::Sadd(npts, -coords[1], ptsy,1,tmpy,1);
662  Vmath::Vmul (npts, tmpx,1,tmpx,1,tmpx,1);
663  Vmath::Vvtvp(npts, tmpy,1,tmpy,1,tmpx,1,tmpx,1);
664 
665  int min_i = Vmath::Imin(npts,tmpx,1);
666 
667  Lcoords[0] = za[min_i%za.num_elements()];
668  Lcoords[1] = zb[min_i/za.num_elements()];
669 
670  // Perform newton iteration to find local coordinates
671  NewtonIterationForLocCoord(coords, ptsx, ptsy, Lcoords,resid);
672  }
673  return resid;
674  }
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:824
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:531
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 786 of file QuadGeom.cpp.

References kNedges.

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 777 of file QuadGeom.cpp.

References kNverts.

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

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 700 of file QuadGeom.cpp.

References ASSERTL2, and m_verts.

701  {
702  ASSERTL2((i >=0) && (i <= 3),"Vertex id must be between 0 and 3");
703  return m_verts[i];
704  }
#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 690 of file QuadGeom.cpp.

References ASSERTL2, and m_verts.

691  {
692  ASSERTL2((i >=0) && (i <= 3),"Verted id must be between 0 and 3");
693  return m_verts[i]->GetVid();
694  }
#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 376 of file QuadGeom.cpp.

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

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 367 of file QuadGeom.cpp.

References m_elmtMap.

368  {
369  return int(m_elmtMap.size());
370  }
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 827 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().

829  {
830  Geometry::v_Reset(curvedEdges, curvedFaces);
831  CurveMap::iterator it = curvedFaces.find(m_globalID);
832 
833  if (it != curvedFaces.end())
834  {
835  m_curve = it->second;
836  }
837 
838  for (int i = 0; i < 4; ++i)
839  {
840  m_edges[i]->Reset(curvedEdges, curvedFaces);
841  }
842 
843  SetUpXmap();
844  SetUpCoeffs(m_xmap->GetNcoeffs());
845  }
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 520 of file QuadGeom.cpp.

References m_ownData.

521  {
522  m_ownData = true;
523  }
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 754 of file QuadGeom.cpp.

References Nektar::iterator, and m_edges.

755  {
756  int returnval = -1;
757 
758  SegGeomVector::iterator edgeIter;
759  int i;
760 
761  for (i=0,edgeIter = m_edges.begin(); edgeIter != m_edges.end(); ++edgeIter,++i)
762  {
763  if (*edgeIter == edge)
764  {
765  returnval = i;
766  break;
767  }
768  }
769 
770  return returnval;
771  }
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.