44 namespace SpatialDomains
47 {0,3,4},{0,1,5},{1,2,6},{2,3,7},
48 {4,8,11},{5,8,9},{6,9,10},{7,10,11}};
50 {0,1,4},{0,1,2},{0,2,3},{0,3,4},
51 {1,4,5},{1,2,5},{2,3,5},{3,4,5}};
53 {0,1},{0,2},{0,3},{0,4},{1,4},{1,2},{2,3},{3,4},
54 {1,5},{2,5},{3,5},{4,5}};
62 Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
98 if (
m_xmap->GetBasisNumModes(0) != 2 ||
99 m_xmap->GetBasisNumModes(1) != 2 ||
100 m_xmap->GetBasisNumModes(2) != 2 )
122 if( fabs( (*
m_verts[ faceVerts[f][0] ])(i) -
123 (*
m_verts[ faceVerts[f][1] ])(i) +
124 (*
m_verts[ faceVerts[f][2] ])(i) -
125 (*
m_verts[ faceVerts[f][3] ])(i) )
171 nq0 =
m_xmap->GetNumPoints(0);
172 nq1 =
m_xmap->GetNumPoints(1);
173 nq2 =
m_xmap->GetNumPoints(2);
178 len0 += (pts[nq0-1]-pts[0])*(pts[nq0-1]-pts[0]);
179 xi0 += (coords[i] -pts[0])*(pts[nq0-1]-pts[0]);
182 len1 += (pts[nq0*(nq1-1)]-pts[0])*(pts[nq0*(nq1-1)]-pts[0]);
183 xi1 += (coords[i] -pts[0])*(pts[nq0*(nq1-1)]-pts[0]);
186 len2 += (pts[nq0*nq1*(nq2-1)]-pts[0])*(pts[nq0*nq1*(nq2-1)]-pts[0]);
187 xi2 += (coords[i] -pts[0])*(pts[nq0*nq1*(nq2-1)]-pts[0]);
190 Lcoords[0] = 2*xi0/len0-1.0;
191 Lcoords[1] = 2*xi1/len1-1.0;
192 Lcoords[2] = 2*xi2/len2-1.0;
196 PointGeom r(m_coordim, 0, coords[0], coords[1], coords[2]);
197 for(
int i = 0; i < 8; ++i)
228 ptdist = sqrt(tmp1[min_i]);
231 int qa = za.num_elements(), qb = zb.num_elements();
232 Lcoords[2] = zc[min_i/(qa*qb)];
233 min_i = min_i%(qa*qb);
234 Lcoords[1] = zb[min_i/qa];
235 Lcoords[0] = za[min_i%qa];
272 ASSERTL1(gloCoord.num_elements() == 3,
273 "Three dimensional geometry expects three coordinates.");
289 for(i = 0; i < 3; ++i)
293 mincoord[i] =
Vmath::Vmin(pts.num_elements(),pts,1);
294 maxcoord[i] =
Vmath::Vmax(pts.num_elements(),pts,1);
296 diff = max(maxcoord[i] - mincoord[i],diff);
299 for(i = 0; i < 3; ++i)
301 if((gloCoord[i] < mincoord[i] - 0.2*diff)||
302 (gloCoord[i] > maxcoord[i] + 0.2*diff))
311 if (locCoord[0] >= -(1+tol) && locCoord[0] <= 1+tol
312 && locCoord[1] >= -(1+tol) && locCoord[1] <= 1+tol
313 && locCoord[2] >= -(1+tol) && locCoord[2] <= 1+tol)
322 for(
int i = 0; i < 3; ++i)
324 if(locCoord[i] <-(1+tol))
326 locCoord[i] = -(1+tol);
329 if(locCoord[i] > (1+tol))
370 if (faceidx == 0 || faceidx == 1)
374 else if (faceidx == 1 || faceidx == 3)
394 for(f = 1; f < 5 ; f++)
397 for(i = 0; i < 4; i++)
399 for(j = 0; j < 4; j++)
412 std::ostringstream errstrm;
413 errstrm <<
"Connected faces do not share an edge. Faces ";
419 std::ostringstream errstrm;
420 errstrm <<
"Connected faces share more than one edge. Faces ";
428 for(i = 0; i < 4; i++)
430 for(j = 0; j < 4; j++)
442 std::ostringstream errstrm;
443 errstrm <<
"Connected faces do not share an edge. Faces ";
449 std::ostringstream errstrm;
450 errstrm <<
"Connected faces share more than one edge. Faces ";
454 for(f = 1; f < 4 ; f++)
457 for(i = 0; i < 4; i++)
459 for(j = 0; j < 4; j++)
472 std::ostringstream errstrm;
473 errstrm <<
"Connected faces do not share an edge. Faces ";
479 std::ostringstream errstrm;
480 errstrm <<
"Connected faces share more than one edge. Faces ";
487 for(f = 1; f < 5 ; f++)
490 for(i = 0; i < 4; i++)
492 for(j = 0; j < 4; j++)
505 std::ostringstream errstrm;
506 errstrm <<
"Connected faces do not share an edge. Faces ";
512 std::ostringstream errstrm;
513 errstrm <<
"Connected faces share more than one edge. Faces ";
537 std::ostringstream errstrm;
538 errstrm <<
"Connected edges do not share a vertex. Edges ";
545 for(i = 1; i < 3; i++)
557 std::ostringstream errstrm;
558 errstrm <<
"Connected edges do not share a vertex. Edges ";
559 errstrm <<
m_edges[i]->GetEid() <<
", " <<
m_edges[i-1]->GetEid();
580 std::ostringstream errstrm;
581 errstrm <<
"Connected edges do not share a vertex. Edges ";
587 for(i = 9; i < 11; i++)
599 std::ostringstream errstrm;
600 errstrm <<
"Connected edges do not share a vertex. Edges ";
601 errstrm <<
m_edges[i]->GetEid() <<
", " <<
m_edges[i-1]->GetEid();
630 unsigned int baseVertex;
653 unsigned int orientation;
659 elementAaxis_length = 0.0;
660 elementBaxis_length = 0.0;
661 faceAaxis_length = 0.0;
662 faceBaxis_length = 0.0;
667 baseVertex =
m_faces[f]->GetVid(0);
679 elementAaxis[i] = (*
m_verts[ faceVerts[f][1] ])[i] - (*
m_verts[ faceVerts[f][0] ])[i];
680 elementBaxis[i] = (*
m_verts[ faceVerts[f][3] ])[i] - (*
m_verts[ faceVerts[f][0] ])[i];
683 else if( baseVertex ==
m_verts[ faceVerts[f][1] ]->
GetVid() )
687 elementAaxis[i] = (*
m_verts[ faceVerts[f][1] ])[i] - (*
m_verts[ faceVerts[f][0] ])[i];
688 elementBaxis[i] = (*
m_verts[ faceVerts[f][2] ])[i] - (*
m_verts[ faceVerts[f][1] ])[i];
691 else if( baseVertex ==
m_verts[ faceVerts[f][2] ]->
GetVid() )
695 elementAaxis[i] = (*
m_verts[ faceVerts[f][2] ])[i] - (*
m_verts[ faceVerts[f][3] ])[i];
696 elementBaxis[i] = (*
m_verts[ faceVerts[f][2] ])[i] - (*
m_verts[ faceVerts[f][1] ])[i];
699 else if( baseVertex ==
m_verts[ faceVerts[f][3] ]->
GetVid() )
703 elementAaxis[i] = (*
m_verts[ faceVerts[f][2] ])[i] - (*
m_verts[ faceVerts[f][3] ])[i];
704 elementBaxis[i] = (*
m_verts[ faceVerts[f][3] ])[i] - (*
m_verts[ faceVerts[f][0] ])[i];
709 ASSERTL0(
false,
"Could not find matching vertex for the face");
719 elementAaxis_length += pow(elementAaxis[i],2);
720 elementBaxis_length += pow(elementBaxis[i],2);
721 faceAaxis_length += pow(faceAaxis[i],2);
722 faceBaxis_length += pow(faceBaxis[i],2);
725 elementAaxis_length = sqrt(elementAaxis_length);
726 elementBaxis_length = sqrt(elementBaxis_length);
727 faceAaxis_length = sqrt(faceAaxis_length);
728 faceBaxis_length = sqrt(faceBaxis_length);
734 dotproduct1 += elementAaxis[i]*faceAaxis[i];
744 if(dotproduct1 < 0.0)
752 dotproduct2 += elementBaxis[i]*faceBaxis[i];
756 ASSERTL1(fabs(elementBaxis_length*faceBaxis_length - fabs(dotproduct2)) <
758 "These vectors should be parallel");
762 if( dotproduct2 < 0.0 )
777 dotproduct1 += elementAaxis[i]*faceBaxis[i];
781 ASSERTL1(fabs(elementAaxis_length*faceBaxis_length - fabs(dotproduct1)) <
783 "These vectors should be parallel");
787 if(dotproduct1 < 0.0)
796 dotproduct2 += elementBaxis[i]*faceAaxis[i];
800 ASSERTL1(fabs(elementBaxis_length*faceAaxis_length - fabs(dotproduct2)) <
802 "These vectors should be parallel");
804 if( dotproduct2 < 0.0 )
810 orientation = orientation + 5;
822 const unsigned int edgeVerts[
kNedges][2] =
849 ASSERTL0(
false,
"Could not find matching vertex for the edge");
860 for (
int i = 0; i < 6; ++i)
862 m_faces[i]->Reset(curvedEdges, curvedFaces);
877 vector<int> tmp1, tmp2;
909 int order0 = *max_element(tmp1.begin(), tmp1.end());
910 int points0 = *max_element(tmp2.begin(), tmp2.end());
945 int order1 = *max_element(tmp1.begin(), tmp1.end());
946 int points1 = *max_element(tmp2.begin(), tmp2.end());
981 int order2 = *max_element(tmp1.begin(), tmp1.end());
982 int points2 = *max_element(tmp2.begin(), tmp2.end());
static const unsigned int VertexEdgeConnectivity[8][3]
StdRegions::StdExpansionSharedPtr m_xmap
void SetUpXmap()
Set up the m_xmap object by determining the order of each direction from derived faces.
#define ASSERTL0(condition, msg)
std::vector< StdRegions::Orientation > m_eorient
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
virtual void v_GenGeomFactors()
GeomFactorsSharedPtr m_geomFactors
T Vmax(int n, const T *x, const int incx)
Return the maximum element in x – called vmax to avoid conflict with max.
static const int kNtfaces
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
boost::shared_ptr< QuadGeom > QuadGeomSharedPtr
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_GetNumFaces() const
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state and remove allocated GeomFactors.
virtual int v_GetNumEdges() const
StdRegions::StdExpansionSharedPtr GetXmap() const
GeomState m_geomFactorsState
std::vector< StdRegions::Orientation > m_forient
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
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...
void SetUpLocalVertices()
boost::shared_ptr< SegGeom > SegGeomSharedPtr
virtual NekDouble v_GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
virtual int v_GetVertexEdgeMap(const int i, const int j) const
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state and remove allocated GeomFactors.
Defines a specification for a set of points.
void Sadd(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Add vector y = alpha + x.
void SetUpEdgeOrientation()
PointGeomSharedPtr GetVertex(int i) const
virtual void v_FillGeom()
Put all quadrature information into face/edge structure and backward transform.
virtual int v_GetEdgeFaceMap(const int i, const int j) const
const Geometry1DSharedPtr GetEdge(int i) const
virtual int v_GetVertexFaceMap(const int i, const int j) const
Array< OneD, Array< OneD, NekDouble > > m_coeffs
static const unsigned int EdgeFaceConnectivity[12][2]
virtual int v_GetNumVerts() const
Geometry is straight-sided with constant geometric factors.
void SetUpFaceOrientation()
virtual int v_GetDir(const int faceidx, const int facedir) const
LibUtilities::ShapeType m_shapeType
static const unsigned int VertexFaceConnectivity[8][3]
NekDouble dist(PointGeom &a)
boost::unordered_map< int, CurveSharedPtr > CurveMap
Geometric information has been generated.
void SetUpCoeffs(const int nCoeffs)
Initialise the m_coeffs array.
static const int kNqfaces
GeomFactorsSharedPtr GetMetricInfo()
GeomType
Indicates the type of element geometry.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Geometry is curved or has non-constant factors.
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.