44 namespace SpatialDomains
47 {0,2,3},{0,1,4},{1,2,5},{3,4,5}};
49 {0,1,3},{0,1,2},{0,2,3},{1,2,3}};
51 {0,1},{0,2},{0,3},{1,3},{1,2},{2,3}};
59 Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
112 ASSERTL1(gloCoord.num_elements() == 3,
113 "Three dimensional geometry expects three coordinates.");
129 for(i = 0; i < 3; ++i)
133 mincoord[i] =
Vmath::Vmin(pts.num_elements(),pts,1);
134 maxcoord[i] =
Vmath::Vmax(pts.num_elements(),pts,1);
136 diff = max(maxcoord[i] - mincoord[i],diff);
139 for(i = 0; i < 3; ++i)
141 if((gloCoord[i] < mincoord[i] - 0.2*diff)||
142 (gloCoord[i] > maxcoord[i] + 0.2*diff))
153 if (locCoord[0] >= -(1+tol) && locCoord[1] >= -(1+tol) &&
154 locCoord[2] >= -(1+tol) &&
155 locCoord[0] + locCoord[1] + locCoord[2] <= -1+tol)
164 for(
int i = 0; i < 3; ++i)
166 if(locCoord[i] <-(1+tol))
168 locCoord[i] = -(1+tol);
171 if(locCoord[i] > (1+tol))
204 cp1020.
Mult(e10,e20);
205 cp2030.
Mult(e20,e30);
206 cp3010.
Mult(e30,e10);
216 Lcoords[0] = 2.0*beta - 1.0;
217 Lcoords[1] = 2.0*gamma - 1.0;
218 Lcoords[2] = 2.0*delta - 1.0;
221 for(
int i = 0; i < 4; ++i)
254 ptdist = sqrt(tmp1[min_i]);
257 int qa = za.num_elements(), qb = zb.num_elements();
258 Lcoords[2] = zc[min_i/(qa*qb)];
259 min_i = min_i%(qa*qb);
260 Lcoords[1] = zb[min_i/qa];
261 Lcoords[0] = za[min_i%qa];
264 Lcoords[1] = (1.0+Lcoords[0])*(1.0-Lcoords[2])/2 -1.0;
265 Lcoords[0] = (1.0+Lcoords[0])*(-Lcoords[1]-Lcoords[2])/2 -1.0;
295 else if (faceidx == 1)
333 std::ostringstream errstrm;
334 errstrm <<
"Local edge 0 (eid=" <<
m_faces[0]->GetEid(0);
335 errstrm <<
") on face " <<
m_faces[0]->GetFid();
336 errstrm <<
" must be the same as local edge 0 (eid="<<
m_faces[1]->GetEid(0);
337 errstrm <<
") on face " <<
m_faces[1]->GetFid();
342 for(faceConnected = 1; faceConnected < 4 ; faceConnected++)
345 for(i = 0; i < 3; i++)
357 std::ostringstream errstrm;
358 errstrm <<
"Face 0 does not share an edge with first edge of adjacent face. Faces ";
364 std::ostringstream errstrm;
365 errstrm <<
"Connected faces share more than one edge. Faces ";
374 for(i = 0; i < 3; i++)
376 for(j = 0; j < 3; j++)
388 std::ostringstream errstrm;
389 errstrm <<
"Connected faces do not share an edge. Faces ";
395 std::ostringstream errstrm;
396 errstrm <<
"Connected faces share more than one edge. Faces ";
401 for(faceConnected = 1; faceConnected < 3 ; faceConnected++)
404 for(i = 0; i < 3; i++)
406 for(j = 0; j < 3; j++)
419 std::ostringstream errstrm;
420 errstrm <<
"Connected faces do not share an edge. Faces ";
426 std::ostringstream errstrm;
427 errstrm <<
"Connected faces share more than one edge. Faces ";
451 std::ostringstream errstrm;
452 errstrm <<
"Connected edges do not share a vertex. Edges ";
458 for(
int i = 1; i < 2; i++)
470 std::ostringstream errstrm;
471 errstrm <<
"Connected edges do not share a vertex. Edges ";
472 errstrm <<
m_edges[i]->GetEid() <<
", " <<
m_edges[i-1]->GetEid();
488 for (
int i = 4; i < 6; ++i)
499 std::ostringstream errstrm;
500 errstrm <<
"Connected edges do not share a vertex. Edges ";
511 const unsigned int edgeVerts[
kNedges][2] =
532 ASSERTL0(
false,
"Could not find matching vertex for the edge");
561 unsigned int baseVertex;
582 unsigned int orientation;
588 elementAaxis_length = 0.0;
589 elementBaxis_length = 0.0;
590 faceAaxis_length = 0.0;
591 faceBaxis_length = 0.0;
596 baseVertex =
m_faces[f]->GetVid(0);
610 elementAaxis[i] = (*
m_verts[ faceVerts[f][1] ])[i] - (*
m_verts[ faceVerts[f][0] ])[i];
611 elementBaxis[i] = (*
m_verts[ faceVerts[f][2] ])[i] - (*
m_verts[ faceVerts[f][0] ])[i];
614 else if( baseVertex ==
m_verts[ faceVerts[f][1] ]->
GetVid() )
618 elementAaxis[i] = (*
m_verts[ faceVerts[f][1] ])[i] - (*
m_verts[ faceVerts[f][0] ])[i];
619 elementBaxis[i] = (*
m_verts[ faceVerts[f][2] ])[i] - (*
m_verts[ faceVerts[f][1] ])[i];
622 else if( baseVertex ==
m_verts[ faceVerts[f][2] ]->
GetVid() )
626 elementAaxis[i] = (*
m_verts[ faceVerts[f][1] ])[i] - (*
m_verts[ faceVerts[f][2] ])[i];
627 elementBaxis[i] = (*
m_verts[ faceVerts[f][2] ])[i] - (*
m_verts[ faceVerts[f][0] ])[i];
632 ASSERTL0(
false,
"Could not find matching vertex for the face");
642 elementAaxis_length += pow(elementAaxis[i],2);
643 elementBaxis_length += pow(elementBaxis[i],2);
644 faceAaxis_length += pow(faceAaxis[i],2);
645 faceBaxis_length += pow(faceBaxis[i],2);
648 elementAaxis_length = sqrt(elementAaxis_length);
649 elementBaxis_length = sqrt(elementBaxis_length);
650 faceAaxis_length = sqrt(faceAaxis_length);
651 faceBaxis_length = sqrt(faceBaxis_length);
657 dotproduct1 += elementAaxis[i]*faceAaxis[i];
667 if(dotproduct1 < 0.0)
675 dotproduct2 += elementBaxis[i]*faceBaxis[i];
679 ASSERTL1(fabs(elementBaxis_length*faceBaxis_length - fabs(dotproduct2)) <
681 "These vectors should be parallel");
685 if( dotproduct2 < 0.0 )
700 dotproduct1 += elementAaxis[i]*faceBaxis[i];
704 ASSERTL1(fabs(elementAaxis_length*faceBaxis_length - fabs(dotproduct1)) <
706 "These vectors should be parallel");
710 if(dotproduct1 < 0.0)
719 dotproduct2 += elementBaxis[i]*faceAaxis[i];
723 ASSERTL1(fabs(elementBaxis_length*faceAaxis_length - fabs(dotproduct2)) <
725 "These vectors should be parallel");
727 if( dotproduct2 < 0.0 )
733 orientation = orientation + 5;
746 for (
int i = 0; i < 4; ++i)
748 m_faces[i]->Reset(curvedEdges, curvedFaces);
763 int order0 = *max_element(tmp.begin(), tmp.end());
767 int points0 = *max_element(tmp.begin(), tmp.end());
770 tmp.push_back(order0);
773 int order1 = *max_element(tmp.begin(), tmp.end());
776 tmp.push_back(points0);
779 int points1 = *max_element(tmp.begin(), tmp.end());
782 tmp.push_back(order0);
783 tmp.push_back(order1);
787 int order2 = *max_element(tmp.begin(), tmp.end());
790 tmp.push_back(points0);
791 tmp.push_back(points1);
795 int points2 = *max_element(tmp.begin(), tmp.end());
virtual int v_GetEdgeFaceMap(const int i, const int j) const
StdRegions::StdExpansionSharedPtr m_xmap
#define ASSERTL0(condition, msg)
void SetUpFaceOrientation()
std::vector< StdRegions::Orientation > m_eorient
Principle Modified Functions .
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
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...
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
NekDouble dot(PointGeom &a)
int Imin(int n, const T *x, const int incx)
Return the index of the minimum element in x.
int GetEid(int i) const
Return the ID of edge i in this element.
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
Principle Modified Functions .
virtual int v_GetDir(const int faceidx, const int facedir) const
void Sub(PointGeom &a, PointGeom &b)
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state and remove allocated GeomFactors.
static const int kNtfaces
StdRegions::StdExpansionSharedPtr GetXmap() const
std::vector< StdRegions::Orientation > m_forient
Gauss Radau pinned at x=-1, .
void SetUpEdgeOrientation()
static const unsigned int EdgeFaceConnectivity[6][2]
virtual int v_GetNumFaces() const
void NewtonIterationForLocCoord(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &ptsx, const Array< OneD, const NekDouble > &ptsy, const Array< OneD, const NekDouble > &ptsz, Array< OneD, NekDouble > &Lcoords, NekDouble &resid)
static const NekDouble kNekZeroTol
boost::shared_ptr< SegGeom > SegGeomSharedPtr
Principle Modified Functions .
virtual int v_GetNumEdges() const
static const unsigned int VertexFaceConnectivity[4][3]
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
Defines a specification for a set of points.
virtual int v_GetVertexEdgeMap(const int i, const int j) const
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
static const int kNqfaces
PointGeomSharedPtr GetVertex(int i) const
virtual void v_FillGeom()
Put all quadrature information into face/edge structure and backward transform.
virtual NekDouble v_GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Get Local cartesian points.
void Mult(PointGeom &a, PointGeom &b)
const Geometry1DSharedPtr GetEdge(int i) const
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state and remove allocated GeomFactors.
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Geometry is straight-sided with constant geometric factors.
virtual int v_GetVertexFaceMap(const int i, const int j) const
Gauss Radau pinned at x=-1, .
virtual int v_GetNumVerts() const
boost::shared_ptr< TriGeom > TriGeomSharedPtr
LibUtilities::ShapeType m_shapeType
NekDouble dist(PointGeom &a)
boost::unordered_map< int, CurveSharedPtr > CurveMap
void SetUpCoeffs(const int nCoeffs)
Initialise the m_coeffs array.
static const unsigned int VertexEdgeConnectivity[4][3]
GeomFactorsSharedPtr GetMetricInfo()
void SetUpLocalVertices()
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
int m_coordim
coordinate dimension
Describes the specification for a Basis.
1D Gauss-Lobatto-Legendre quadrature points
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.