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

#include <TriGeom.h>

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

Public Member Functions

 TriGeom ()
 
 TriGeom (int id, const int coordim)
 
 TriGeom (const int id, const PointGeomSharedPtr verts[], const SegGeomSharedPtr edges[], const StdRegions::Orientation eorient[])
 
 TriGeom (const int id, const SegGeomSharedPtr edges[], const StdRegions::Orientation eorient[])
 
 TriGeom (const int id, const SegGeomSharedPtr edges[], const StdRegions::Orientation eorient[], const CurveSharedPtr &curve)
 
 TriGeom (const TriGeom &in)
 
 ~TriGeom ()
 
NekDouble GetCoord (const int i, const Array< OneD, const NekDouble > &Lcoord)
 
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry2D
 Geometry2D ()
 
 Geometry2D (const int coordim)
 
virtual ~Geometry2D ()
 
int GetFid () const
 
const Geometry1DSharedPtr GetEdge (int i) const
 
const Geometry2DSharedPtr GetFace (int i) const
 
int WhichEdge (SegGeomSharedPtr edge)
 
int WhichFace (Geometry2DSharedPtr face)
 
const LibUtilities::BasisSharedPtr GetEdgeBasis (const int i)
 
StdRegions::Orientation GetFaceOrient (const int i) const
 
StdRegions::Orientation GetCartesianEorient (const int i) const
 
- Public Member Functions inherited from Nektar::SpatialDomains::Geometry
 Geometry ()
 
 Geometry (int coordim)
 
virtual ~Geometry ()
 
bool IsElmtConnected (int gvo_id, int locid) const
 
void AddElmtConnected (int gvo_id, int locid)
 
int NumElmtConnected () const
 
int GetCoordim () const
 
void SetCoordim (int coordim)
 
GeomFactorsSharedPtr GetGeomFactors ()
 
GeomFactorsSharedPtr GetRefGeomFactors (const Array< OneD, const LibUtilities::BasisSharedPtr > &tbasis)
 
GeomFactorsSharedPtr GetMetricInfo ()
 
LibUtilities::ShapeType GetShapeType (void)
 
int GetGlobalID (void)
 
void SetGlobalID (int globalid)
 
int GetVid (int i) const
 
int GetEid (int i) const
 
int GetFid (int i) const
 
int GetTid (int i) const
 
int GetNumVerts () const
 
PointGeomSharedPtr GetVertex (int i) const
 
StdRegions::Orientation GetEorient (const int i) const
 
StdRegions::Orientation GetPorient (const int i) const
 
StdRegions::Orientation GetForient (const int i) const
 
int GetNumEdges () const
 
int GetNumFaces () const
 
int GetShapeDim () const
 
StdRegions::StdExpansionSharedPtr GetXmap () const
 
const Array< OneD, const
NekDouble > & 
GetCoeffs (const int i) const
 
bool ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
 
bool ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol)
 
bool ContainsPoint (const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &resid)
 
int GetVertexEdgeMap (int i, int j) const
 
int GetVertexFaceMap (int i, int j) const
 return the id of the $j^{th}$ face attached to the $ i^{th}$ vertex 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 TriGeom &face1, const TriGeom &face2)
 
static StdRegions::Orientation GetFaceOrientation (const PointGeomVector &face1, const PointGeomVector &face2)
 

Static Public Attributes

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

Protected Member Functions

virtual void v_AddElmtConnected (int gvo_id, int locid)
 
virtual int v_NumElmtConnected () const
 
virtual bool v_IsElmtConnected (int gvo_id, int locid) const
 
virtual int v_GetFid () const
 
virtual int v_GetCoordim () const
 
virtual const
LibUtilities::BasisSharedPtr 
v_GetBasis (const int i)
 
virtual const
LibUtilities::BasisSharedPtr 
v_GetEdgeBasis (const int i)
 
virtual NekDouble v_GetCoord (const int i, const Array< OneD, const NekDouble > &Lcoord)
 
virtual void v_GenGeomFactors ()
 
virtual void v_SetOwnData ()
 
virtual void v_FillGeom ()
 Put all quadrature information into edge structure. 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)
 Determines if a point specified in global coordinates is located within this tetrahedral geometry. More...
 
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
 

Additional Inherited Members

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

Detailed Description

Definition at line 65 of file TriGeom.h.

Constructor & Destructor Documentation

Nektar::SpatialDomains::TriGeom::TriGeom ( )

Definition at line 54 of file TriGeom.cpp.

References Nektar::LibUtilities::eTriangle.

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

Definition at line 62 of file TriGeom.cpp.

References m_fid, Nektar::SpatialDomains::Geometry::m_globalID, Nektar::SpatialDomains::Geometry::m_xmap, Nektar::SpatialDomains::Geometry::SetUpCoeffs(), and SetUpXmap().

62  :
63  Geometry2D(coordim), m_fid(id)
64  {
65  m_globalID = m_fid;
66  SetUpXmap();
67  SetUpCoeffs(m_xmap->GetNcoeffs());
68  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
void SetUpCoeffs(const int nCoeffs)
Initialise the m_coeffs array.
Definition: Geometry.h:484
Nektar::SpatialDomains::TriGeom::TriGeom ( const int  id,
const PointGeomSharedPtr  verts[],
const SegGeomSharedPtr  edges[],
const StdRegions::Orientation  eorient[] 
)

Copy the vert shared pointers.

Copy the edge shared pointers.

Definition at line 73 of file TriGeom.cpp.

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

76  :
77  Geometry2D(verts[0]->GetCoordim()),
78  m_fid(id)
79  {
80  m_globalID = m_fid;
82 
83  /// Copy the vert shared pointers.
84  m_verts.insert(m_verts.begin(), verts, verts+TriGeom::kNverts);
85 
86  /// Copy the edge shared pointers.
87  m_edges.insert(m_edges.begin(), edges, edges+TriGeom::kNedges);
88 
89  for (int j=0; j<kNedges; ++j)
90  {
91  m_eorient[j] = eorient[j];
92  }
93 
94  m_coordim = verts[0]->GetCoordim();
95  ASSERTL0(m_coordim > 1,
96  "Cannot call function with dim == 1");
97 
98  SetUpXmap();
99  SetUpCoeffs(m_xmap->GetNcoeffs());
100  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
static const int kNverts
Definition: TriGeom.h:99
StdRegions::Orientation m_eorient[kNedges]
Definition: TriGeom.h:113
static const int kNedges
Get the orientation of face1.
Definition: TriGeom.h:98
LibUtilities::ShapeType m_shapeType
Definition: Geometry.h:177
void SetUpCoeffs(const int nCoeffs)
Initialise the m_coeffs array.
Definition: Geometry.h:484
int m_coordim
coordinate dimension
Definition: Geometry.h:169
Nektar::SpatialDomains::TriGeom::TriGeom ( const int  id,
const SegGeomSharedPtr  edges[],
const StdRegions::Orientation  eorient[] 
)

Copy the edge shared pointers.

Definition at line 106 of file TriGeom.cpp.

References ASSERTL0, Nektar::StdRegions::eForwards, Nektar::LibUtilities::eTriangle, 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().

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

Copy the edge shared pointers.

Definition at line 145 of file TriGeom.cpp.

References ASSERTL0, Nektar::StdRegions::eForwards, Nektar::LibUtilities::eTriangle, 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().

150  : Geometry2D(edges[0]->GetVertex(0)->GetCoordim()),
151  m_fid(id),
152  m_curve(curve)
153  {
154  m_globalID = m_fid;
156 
157  /// Copy the edge shared pointers.
158  m_edges.insert(m_edges.begin(), edges, edges+TriGeom::kNedges);
159 
160  for(int j=0; j <kNedges; ++j)
161  {
162  if(eorient[j] == StdRegions::eForwards)
163  {
164  m_verts.push_back(edges[j]->GetVertex(0));
165  }
166  else
167  {
168  m_verts.push_back(edges[j]->GetVertex(1));
169  }
170  }
171 
172  for (int j=0; j<kNedges; ++j)
173  {
174  m_eorient[j] = eorient[j];
175  }
176 
177  m_coordim = edges[0]->GetVertex(0)->GetCoordim();
178  ASSERTL0(m_coordim > 1,"Cannot call function with dim == 1");
179 
180  SetUpXmap();
181  SetUpCoeffs(m_xmap->GetNcoeffs());
182  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
StdRegions::Orientation m_eorient[kNedges]
Definition: TriGeom.h:113
PointGeomSharedPtr GetVertex(int i) const
Definition: Geometry.h:348
static const int kNedges
Get the orientation of face1.
Definition: TriGeom.h:98
LibUtilities::ShapeType m_shapeType
Definition: Geometry.h:177
void SetUpCoeffs(const int nCoeffs)
Initialise the m_coeffs array.
Definition: Geometry.h:484
int m_coordim
coordinate dimension
Definition: Geometry.h:169
Nektar::SpatialDomains::TriGeom::TriGeom ( const TriGeom in)

Definition at line 187 of file TriGeom.cpp.

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

188  {
189  // From Geomtry
190  m_shapeType = in.m_shapeType;
191 
192  // From TriFaceComponent
193  m_fid = in.m_fid;
194  m_globalID = m_fid;
195  m_ownVerts = in.m_ownVerts;
196  std::list<CompToElmt>::const_iterator def;
197  for(def = in.m_elmtMap.begin(); def != in.m_elmtMap.end(); def++)
198  {
199  m_elmtMap.push_back(*def);
200  }
201 
202  // From TriGeom
203  m_verts = in.m_verts;
204  m_edges = in.m_edges;
205  for (int i = 0; i < kNedges; i++)
206  {
207  m_eorient[i] = in.m_eorient[i];
208  }
209  m_ownData = in.m_ownData;
210  }
StdRegions::Orientation m_eorient[kNedges]
Definition: TriGeom.h:113
std::list< CompToElmt > m_elmtMap
Definition: TriGeom.h:116
static const int kNedges
Get the orientation of face1.
Definition: TriGeom.h:98
LibUtilities::ShapeType m_shapeType
Definition: Geometry.h:177
Nektar::SpatialDomains::TriGeom::~TriGeom ( )

Definition at line 216 of file TriGeom.cpp.

217  {
218  }

Member Function Documentation

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

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

Definition at line 225 of file TriGeom.cpp.

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

Referenced by v_GetCoord().

227  {
229  "Geometry is not in physical space");
230 
231  Array<OneD, NekDouble> tmp(m_xmap->GetTotPoints());
232  m_xmap->BwdTrans(m_coeffs[i], tmp);
233 
234  return m_xmap->PhysEvaluate(Lcoord, tmp);
235  }
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:218
StdRegions::Orientation Nektar::SpatialDomains::TriGeom::GetFaceOrientation ( const TriGeom face1,
const TriGeom face2 
)
static

Definition at line 237 of file TriGeom.cpp.

References m_verts.

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

240  {
241  return GetFaceOrientation(face1.m_verts, face2.m_verts);
242  }
static StdRegions::Orientation GetFaceOrientation(const TriGeom &face1, const TriGeom &face2)
Definition: TriGeom.cpp:237
StdRegions::Orientation Nektar::SpatialDomains::TriGeom::GetFaceOrientation ( const PointGeomVector face1,
const PointGeomVector face2 
)
static

Definition at line 244 of file TriGeom.cpp.

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

247  {
248  int i, j, vmap[3] = {-1,-1,-1};
249  NekDouble x, y, z, x1, y1, z1, cx = 0.0, cy = 0.0, cz = 0.0;
250 
251  // For periodic faces, we calculate the vector between the centre
252  // points of the two faces. (For connected faces this will be
253  // zero). We can then use this to determine alignment later in the
254  // algorithm.
255  for (i = 0; i < 3; ++i)
256  {
257  cx += (*face2[i])(0) - (*face1[i])(0);
258  cy += (*face2[i])(1) - (*face1[i])(1);
259  cz += (*face2[i])(2) - (*face1[i])(2);
260  }
261  cx /= 3;
262  cy /= 3;
263  cz /= 3;
264 
265  // Now construct a mapping which takes us from the vertices of one
266  // face to the other. That is, vertex j of face2 corresponds to
267  // vertex vmap[j] of face1.
268  for (i = 0; i < 3; ++i)
269  {
270  x = (*face1[i])(0);
271  y = (*face1[i])(1);
272  z = (*face1[i])(2);
273  for (j = 0; j < 3; ++j)
274  {
275  x1 = (*face2[j])(0)-cx;
276  y1 = (*face2[j])(1)-cy;
277  z1 = (*face2[j])(2)-cz;
278  if (sqrt((x1-x)*(x1-x)+(y1-y)*(y1-y)+(z1-z)*(z1-z)) < 1e-8)
279  {
280  vmap[j] = i;
281  break;
282  }
283  }
284  }
285 
286  if (vmap[1] == (vmap[0]+1) % 3)
287  {
288  switch (vmap[0])
289  {
290  case 0:
292  break;
293  case 1:
295  break;
296  case 2:
298  break;
299  }
300  }
301  else
302  {
303  switch (vmap[0])
304  {
305  case 0:
307  break;
308  case 1:
310  break;
311  case 2:
313  break;
314  }
315  }
316 
317  ASSERTL0(false, "Unable to determine triangle orientation");
319  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
double NekDouble
void Nektar::SpatialDomains::TriGeom::SetUpXmap ( )
private

Definition at line 850 of file TriGeom.cpp.

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

Referenced by TriGeom(), and v_Reset().

851  {
852  int order0 = m_edges[0]->GetBasis(0)->GetNumModes();
853  int order1 = max(order0,
854  max(m_edges[1]->GetBasis(0)->GetNumModes(),
855  m_edges[2]->GetBasis(0)->GetNumModes()));
856 
857  const LibUtilities::BasisKey B0(
858  LibUtilities::eModified_A, order0, LibUtilities::PointsKey(
860  const LibUtilities::BasisKey B1(
861  LibUtilities::eModified_B, order1, LibUtilities::PointsKey(
863 
865  ::AllocateSharedPtr(B0,B1);
866  }
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
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:57
const LibUtilities::BasisSharedPtr GetBasis(const int i)
Return the j-th basis of the i-th co-ordinate dimension.
Definition: Geometry.h:475
Principle Modified Functions .
Definition: BasisType.h:50
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50
void Nektar::SpatialDomains::TriGeom::v_AddElmtConnected ( int  gvo_id,
int  locid 
)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 325 of file TriGeom.cpp.

References m_elmtMap.

326  {
327  CompToElmt ee(gvo_id,locid);
328  m_elmtMap.push_back(ee);
329  }
std::list< CompToElmt > m_elmtMap
Definition: TriGeom.h:116
bool Nektar::SpatialDomains::TriGeom::v_ContainsPoint ( const Array< OneD, const NekDouble > &  gloCoord,
NekDouble  tol = 0.0 
)
protectedvirtual

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

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 793 of file TriGeom.cpp.

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

Referenced by v_ContainsPoint().

795  {
796  Array<OneD,NekDouble> locCoord(GetCoordim(),0.0);
797  return v_ContainsPoint(gloCoord,locCoord,tol);
798 
799  }
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
Determines if a point specified in global coordinates is located within this tetrahedral geometry...
Definition: TriGeom.cpp:793
bool Nektar::SpatialDomains::TriGeom::v_ContainsPoint ( const Array< OneD, const NekDouble > &  gloCoord,
Array< OneD, NekDouble > &  locCoord,
NekDouble  tol 
)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 804 of file TriGeom.cpp.

References v_ContainsPoint().

807  {
808  NekDouble resid;
809  return v_ContainsPoint(gloCoord,stdCoord,tol,resid);
810  }
double NekDouble
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0)
Determines if a point specified in global coordinates is located within this tetrahedral geometry...
Definition: TriGeom.cpp:793
bool Nektar::SpatialDomains::TriGeom::v_ContainsPoint ( const Array< OneD, const NekDouble > &  gloCoord,
Array< OneD, NekDouble > &  locCoord,
NekDouble  tol,
NekDouble resid 
)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 812 of file TriGeom.cpp.

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

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

Put all quadrature information into edge structure.

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 452 of file TriGeom.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, 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, npts, and Nektar::LibUtilities::PointsManager().

Referenced by v_GenGeomFactors(), and v_GetLocCoords().

453  {
454  // check to see if geometry structure is already filled
455  if(m_state != ePtsFilled)
456  {
457  int i,j,k;
458  int nEdgeCoeffs = m_xmap->GetEdgeNcoeffs(0);
459 
460  if (m_curve)
461  {
462  int pdim = LibUtilities::PointsManager()[
463  LibUtilities::PointsKey(2, m_curve->m_ptype)]
464  ->GetPointsDim();
465 
466  // Deal with 2D points type separately
467  // (e.g. electrostatic or Fekete points) to 1D tensor
468  // product.
469  if (pdim == 2)
470  {
471  int N = m_curve->m_points.size();
472  int nEdgePts = (
473  -1+(int)sqrt(static_cast<NekDouble>(8*N+1)))/2;
474 
475  ASSERTL0(nEdgePts*(nEdgePts+1)/2 == N,
476  "NUMPOINTS should be a triangle number for"
477  " triangle "
478  + boost::lexical_cast<string>(m_globalID));
479 
480  for (i = 0; i < kNedges; ++i)
481  {
482  ASSERTL0(
483  m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
484  "Number of edge points does not correspond "
485  "to number of face points in triangle "
486  + boost::lexical_cast<string>(m_globalID));
487  }
488 
489  const LibUtilities::PointsKey P0(
491  const LibUtilities::PointsKey P1(
493  const LibUtilities::BasisKey T0(
494  LibUtilities::eOrtho_A, nEdgePts, P0);
495  const LibUtilities::BasisKey T1(
496  LibUtilities::eOrtho_B, nEdgePts, P1);
497  Array<OneD, NekDouble> phys(
498  max(nEdgePts*nEdgePts, m_xmap->GetTotPoints()));
499  Array<OneD, NekDouble> tmp(nEdgePts*nEdgePts);
500 
501  for(i = 0; i < m_coordim; ++i)
502  {
503  // Create a StdNodalTriExp.
506  ::AllocateSharedPtr(T0, T1, m_curve->m_ptype);
507 
508  for (j = 0; j < N; ++j)
509  {
510  phys[j] = (m_curve->m_points[j]->GetPtr())[i];
511  }
512 
513  t->BwdTrans(phys, tmp);
514 
515  // Interpolate points to standard region.
517  P0, P1, tmp,
518  m_xmap->GetBasis(0)->GetPointsKey(),
519  m_xmap->GetBasis(1)->GetPointsKey(),
520  phys);
521 
522  // Forwards transform to get coefficient space.
523  m_xmap->FwdTrans(phys, m_coeffs[i]);
524  }
525  }
526  else if (pdim == 1)
527  {
528  int npts = m_curve->m_points.size();
529  int nEdgePts = (int)sqrt(static_cast<NekDouble>(npts));
530  Array<OneD, NekDouble> tmp (npts);
531  Array<OneD, NekDouble> phys(m_xmap->GetTotPoints());
532  LibUtilities::PointsKey curveKey(
533  nEdgePts, m_curve->m_ptype);
534 
535  // Sanity checks:
536  // - Curved faces should have square number of points;
537  // - Each edge should have sqrt(npts) points.
538  ASSERTL0(nEdgePts * nEdgePts == npts,
539  "NUMPOINTS should be a square number for"
540  " triangle "
541  + boost::lexical_cast<string>(m_globalID));
542 
543  for (i = 0; i < kNedges; ++i)
544  {
545  ASSERTL0(
546  m_edges[i]->GetXmap()->GetNcoeffs() == nEdgePts,
547  "Number of edge points does not correspond to "
548  "number of face points in triangle "
549  + boost::lexical_cast<string>(m_globalID));
550  }
551 
552  for (i = 0; i < m_coordim; ++i)
553  {
554  for (j = 0; j < npts; ++j)
555  {
556  tmp[j] = (m_curve->m_points[j]->GetPtr())[i];
557  }
558 
559  // Interpolate curve points to standard triangle
560  // points.
562  curveKey, curveKey, tmp,
563  m_xmap->GetBasis(0)->GetPointsKey(),
564  m_xmap->GetBasis(1)->GetPointsKey(),
565  phys);
566 
567  // Forwards transform to get coefficient space.
568  m_xmap->FwdTrans(phys, m_coeffs[i]);
569  }
570  }
571  else
572  {
573  ASSERTL0(false, "Only 1D/2D points distributions "
574  "supported.");
575  }
576  }
577 
578  Array<OneD, unsigned int> mapArray (nEdgeCoeffs);
579  Array<OneD, int> signArray(nEdgeCoeffs);
580 
581  for(i = 0; i < kNedges; i++)
582  {
583  m_edges[i]->FillGeom();
584  m_xmap->GetEdgeToElementMap(i,m_eorient[i],mapArray,signArray);
585 
586  nEdgeCoeffs = m_edges[i]->GetXmap()->GetNcoeffs();
587 
588  for(j = 0 ; j < m_coordim; j++)
589  {
590  for(k = 0; k < nEdgeCoeffs; k++)
591  {
592  m_coeffs[j][mapArray[k]] =
593  signArray[k] * m_edges[i]->GetCoeffs(j)[k];
594  }
595  }
596  }
597 
599  }
600  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:188
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
StdRegions::StdExpansionSharedPtr GetXmap() const
Definition: Geometry.h:383
boost::shared_ptr< StdNodalTriExp > StdNodalTriExpSharedPtr
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:57
StdRegions::Orientation m_eorient[kNedges]
Definition: TriGeom.h:113
Principle Orthogonal Functions .
Definition: BasisType.h:47
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
PointsManagerT & PointsManager(void)
Principle Orthogonal Functions .
Definition: BasisType.h:46
static const int kNedges
Get the orientation of face1.
Definition: TriGeom.h:98
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
int m_coordim
coordinate dimension
Definition: Geometry.h:169
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50
void Nektar::SpatialDomains::TriGeom::v_GenGeomFactors ( )
protectedvirtual

Set up GeoFac for this geometry using Coord quadrature distribution

Implements Nektar::SpatialDomains::Geometry.

Definition at line 413 of file TriGeom.cpp.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::SpatialDomains::eDeformed, Nektar::SpatialDomains::ePtsFilled, Nektar::SpatialDomains::eRegular, Nektar::SpatialDomains::Geometry::m_coeffs, Nektar::SpatialDomains::Geometry::m_coordim, Nektar::SpatialDomains::Geometry::m_geomFactors, Nektar::SpatialDomains::Geometry::m_geomFactorsState, Nektar::SpatialDomains::Geometry::m_xmap, and v_FillGeom().

414  {
416  {
417  GeomType Gtype = eRegular;
418 
420 
421  // check to see if expansions are linear
422  for(int i = 0; i < m_coordim; ++i)
423  {
424  if(m_xmap->GetBasisNumModes(0) != 2 ||
425  m_xmap->GetBasisNumModes(1) != 2)
426  {
427  Gtype = eDeformed;
428  }
429  }
430 
432  Gtype, m_coordim, m_xmap, m_coeffs);
433 
435  }
436  }
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
virtual void v_FillGeom()
Put all quadrature information into edge structure.
Definition: TriGeom.cpp:452
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::TriGeom::v_GetBasis ( const int  i)
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 377 of file TriGeom.cpp.

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

378  {
379  return m_xmap->GetBasis(i);
380  }
StdRegions::StdExpansionSharedPtr m_xmap
Definition: Geometry.h:172
StdRegions::Orientation Nektar::SpatialDomains::TriGeom::v_GetCartesianEorient ( const int  i) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 727 of file TriGeom.cpp.

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

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 404 of file TriGeom.cpp.

References GetCoord().

405  {
406  return GetCoord(i,Lcoord);
407  }
NekDouble GetCoord(const int i, const Array< OneD, const NekDouble > &Lcoord)
Definition: TriGeom.cpp:225
int Nektar::SpatialDomains::TriGeom::v_GetCoordim ( void  ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 368 of file TriGeom.cpp.

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

369  {
370  return m_coordim;
371  }
int m_coordim
coordinate dimension
Definition: Geometry.h:169
const Geometry1DSharedPtr Nektar::SpatialDomains::TriGeom::v_GetEdge ( int  i) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 707 of file TriGeom.cpp.

References ASSERTL2, and m_edges.

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

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 386 of file TriGeom.cpp.

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

387  {
388  ASSERTL1(i <= 2, "edge is out of range");
389 
390  if(i == 0)
391  {
392  return m_xmap->GetBasis(0);
393  }
394  else
395  {
396  return m_xmap->GetBasis(1);
397  }
398  }
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:218
int Nektar::SpatialDomains::TriGeom::v_GetEid ( int  i) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 677 of file TriGeom.cpp.

References ASSERTL2, and m_edges.

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

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 717 of file TriGeom.cpp.

References ASSERTL2, and m_eorient.

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

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 359 of file TriGeom.cpp.

References m_fid.

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 606 of file TriGeom.cpp.

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

608  {
609  NekDouble resid = 0.0;
611 
612  // calculate local coordinate for coord
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-2 edges.
620  dv1.Sub(*m_verts[1],*m_verts[0]);
621  dv2.Sub(*m_verts[2],*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 proportional to normal of 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  {
642  // Determine nearest point of coords to values in m_xmap
643  int npts = m_xmap->GetTotPoints();
644  Array<OneD, NekDouble> ptsx(npts), ptsy(npts);
645  Array<OneD, NekDouble> tmpx(npts), tmpy(npts);
646 
647  m_xmap->BwdTrans(m_coeffs[0], ptsx);
648  m_xmap->BwdTrans(m_coeffs[1], ptsy);
649 
650  const Array<OneD, const NekDouble> za = m_xmap->GetPoints(0);
651  const Array<OneD, const NekDouble> zb = m_xmap->GetPoints(1);
652 
653  //guess the first local coords based on nearest point
654  Vmath::Sadd(npts, -coords[0], ptsx,1,tmpx,1);
655  Vmath::Sadd(npts, -coords[1], ptsy,1,tmpy,1);
656  Vmath::Vmul (npts, tmpx,1,tmpx,1,tmpx,1);
657  Vmath::Vvtvp(npts, tmpy,1,tmpy,1,tmpx,1,tmpx,1);
658 
659  int min_i = Vmath::Imin(npts,tmpx,1);
660 
661  Lcoords[0] = za[min_i%za.num_elements()];
662  Lcoords[1] = zb[min_i/za.num_elements()];
663 
664  // recover cartesian coordinate from collapsed coordinate.
665  Lcoords[0] = (1.0+Lcoords[0])*(1.0-Lcoords[1])/2 -1.0;
666 
667  // Perform newton iteration to find local coordinates
668  NewtonIterationForLocCoord(coords, ptsx, ptsy, Lcoords,resid);
669  }
670  return resid;
671  }
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
virtual void v_FillGeom()
Put all quadrature information into edge structure.
Definition: TriGeom.cpp:452
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:103
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
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::TriGeom::v_GetNumEdges ( ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 783 of file TriGeom.cpp.

References kNedges.

784  {
785  return kNedges;
786  }
static const int kNedges
Get the orientation of face1.
Definition: TriGeom.h:98
int Nektar::SpatialDomains::TriGeom::v_GetNumVerts ( ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 774 of file TriGeom.cpp.

References kNverts.

775  {
776  return kNverts;
777  }
static const int kNverts
Definition: TriGeom.h:99
PointGeomSharedPtr Nektar::SpatialDomains::TriGeom::v_GetVertex ( int  i) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 697 of file TriGeom.cpp.

References ASSERTL2, and m_verts.

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 687 of file TriGeom.cpp.

References ASSERTL2, and m_verts.

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

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 344 of file TriGeom.cpp.

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

345  {
346  std::list<CompToElmt>::const_iterator def;
347  CompToElmt ee(gvo_id,locid);
348 
349  def = find(m_elmtMap.begin(),m_elmtMap.end(),ee);
350 
351  // Found the element connectivity object in the list
352  return (def != m_elmtMap.end());
353  }
std::list< CompToElmt > m_elmtMap
Definition: TriGeom.h:116
InputIterator find(InputIterator first, InputIterator last, InputIterator startingpoint, const EqualityComparable &value)
Definition: StdRegions.hpp:315
int Nektar::SpatialDomains::TriGeom::v_NumElmtConnected ( ) const
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 335 of file TriGeom.cpp.

References m_elmtMap.

336  {
337  return int(m_elmtMap.size());
338  }
std::list< CompToElmt > m_elmtMap
Definition: TriGeom.h:116
void Nektar::SpatialDomains::TriGeom::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 829 of file TriGeom.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().

832  {
833  Geometry::v_Reset(curvedEdges, curvedFaces);
834  CurveMap::iterator it = curvedFaces.find(m_globalID);
835 
836  if (it != curvedFaces.end())
837  {
838  m_curve = it->second;
839  }
840 
841  for (int i = 0; i < 3; ++i)
842  {
843  m_edges[i]->Reset(curvedEdges, curvedFaces);
844  }
845 
846  SetUpXmap();
847  SetUpCoeffs(m_xmap->GetNcoeffs());
848  }
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::TriGeom::v_SetOwnData ( )
protectedvirtual

Reimplemented from Nektar::SpatialDomains::Geometry.

Definition at line 441 of file TriGeom.cpp.

References m_ownData.

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

Return the edge number of the given edge.

Reimplemented from Nektar::SpatialDomains::Geometry2D.

Definition at line 751 of file TriGeom.cpp.

References Nektar::iterator, and m_edges.

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

Member Data Documentation

const int Nektar::SpatialDomains::TriGeom::kNedges = 3
static
const int Nektar::SpatialDomains::TriGeom::kNverts = 3
static
CurveSharedPtr Nektar::SpatialDomains::TriGeom::m_curve
protected

Definition at line 117 of file TriGeom.h.

Referenced by v_FillGeom(), and v_Reset().

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

Definition at line 112 of file TriGeom.h.

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

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

Definition at line 116 of file TriGeom.h.

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

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

Definition at line 113 of file TriGeom.h.

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

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

Definition at line 114 of file TriGeom.h.

Referenced by TriGeom(), and v_GetFid().

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

Definition at line 198 of file TriGeom.h.

Referenced by TriGeom(), and v_SetOwnData().

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

Definition at line 115 of file TriGeom.h.

Referenced by TriGeom().

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

Definition at line 111 of file TriGeom.h.

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