45 namespace SpatialDomains
47 const unsigned int TetGeom::VertexEdgeConnectivity[4][3] = {
48 {0, 2, 3}, {0, 1, 4}, {1, 2, 5}, {3, 4, 5}};
49 const unsigned int TetGeom::VertexFaceConnectivity[4][3] = {
50 {0, 1, 3}, {0, 1, 2}, {0, 2, 3}, {1, 2, 3}};
51 const unsigned int TetGeom::EdgeFaceConnectivity[6][2] = {
52 {0, 1}, {0, 2}, {0, 3}, {1, 3}, {1, 2}, {2, 3}};
60 :
Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
105 if (locCoord[0] >= -(1 + tol) && locCoord[1] >= -(1 + tol) &&
106 locCoord[2] >= -(1 + tol) &&
107 locCoord[0] + locCoord[1] + locCoord[2] <= -1 + tol)
138 cp1020.
Mult(e10, e20);
139 cp2030.
Mult(e20, e30);
140 cp3010.
Mult(e30, e10);
149 Lcoords[0] = 2.0 * beta - 1.0;
150 Lcoords[1] = 2.0 * gamma - 1.0;
151 Lcoords[2] = 2.0 * delta - 1.0;
154 for (
int i = 0; i < 4; ++i)
164 int npts =
m_xmap->GetTotPoints();
187 ptdist = sqrt(tmp1[min_i]);
190 int qa = za.num_elements(), qb = zb.num_elements();
191 Lcoords[2] = zc[min_i / (qa * qb)];
192 min_i = min_i % (qa * qb);
193 Lcoords[1] = zb[min_i / qa];
194 Lcoords[0] = za[min_i % qa];
197 Lcoords[1] = (1.0 + Lcoords[0]) * (1.0 - Lcoords[2]) / 2 - 1.0;
198 Lcoords[0] = (1.0 + Lcoords[0]) * (-Lcoords[1] - Lcoords[2]) / 2 - 1.0;
213 else if (faceidx == 1)
251 std::ostringstream errstrm;
252 errstrm <<
"Local edge 0 (eid=" <<
m_faces[0]->GetEid(0);
253 errstrm <<
") on face " <<
m_faces[0]->GetGlobalID();
254 errstrm <<
" must be the same as local edge 0 (eid=" 256 errstrm <<
") on face " <<
m_faces[1]->GetGlobalID();
261 for (faceConnected = 1; faceConnected < 4; faceConnected++)
264 for (i = 0; i < 3; i++)
268 edge = dynamic_pointer_cast<
SegGeom>(
277 std::ostringstream errstrm;
278 errstrm <<
"Face 0 does not share an edge with first edge of " 279 "adjacent face. Faces ";
286 std::ostringstream errstrm;
287 errstrm <<
"Connected faces share more than one edge. Faces ";
296 for (i = 0; i < 3; i++)
298 for (j = 0; j < 3; j++)
302 edge = dynamic_pointer_cast<
SegGeom>(
311 std::ostringstream errstrm;
312 errstrm <<
"Connected faces do not share an edge. Faces ";
319 std::ostringstream errstrm;
320 errstrm <<
"Connected faces share more than one edge. Faces ";
326 for (faceConnected = 1; faceConnected < 3; faceConnected++)
329 for (i = 0; i < 3; i++)
331 for (j = 0; j < 3; j++)
333 if ((
m_faces[faceConnected])->GetEid(i) ==
336 edge = dynamic_pointer_cast<
SegGeom>(
346 std::ostringstream errstrm;
347 errstrm <<
"Connected faces do not share an edge. Faces ";
354 std::ostringstream errstrm;
355 errstrm <<
"Connected faces share more than one edge. Faces ";
381 std::ostringstream errstrm;
382 errstrm <<
"Connected edges do not share a vertex. Edges ";
383 errstrm <<
m_edges[0]->GetGlobalID() <<
", " 389 for (
int i = 1; i < 2; i++)
401 std::ostringstream errstrm;
402 errstrm <<
"Connected edges do not share a vertex. Edges ";
403 errstrm <<
m_edges[i]->GetGlobalID() <<
", " 404 <<
m_edges[i - 1]->GetGlobalID();
421 for (
int i = 4; i < 6; ++i)
433 std::ostringstream errstrm;
434 errstrm <<
"Connected edges do not share a vertex. Edges ";
435 errstrm <<
m_edges[3]->GetGlobalID() <<
", " 447 const unsigned int edgeVerts[
kNedges][2] = {
448 {0, 1}, {1, 2}, {0, 2}, {0, 3}, {1, 3}, {2, 3}};
463 ASSERTL0(
false,
"Could not find matching vertex for the edge");
492 unsigned int baseVertex;
505 {0, 1, 2}, {0, 1, 3}, {1, 2, 3}, {0, 2, 3}};
510 unsigned int orientation;
516 elementAaxis_length = 0.0;
517 elementBaxis_length = 0.0;
518 faceAaxis_length = 0.0;
519 faceBaxis_length = 0.0;
524 baseVertex =
m_faces[f]->GetVid(0);
542 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
543 (*
m_verts[faceVerts[f][0]])[i];
544 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
545 (*
m_verts[faceVerts[f][0]])[i];
552 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
553 (*
m_verts[faceVerts[f][0]])[i];
554 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
555 (*
m_verts[faceVerts[f][1]])[i];
562 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
563 (*
m_verts[faceVerts[f][2]])[i];
564 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
565 (*
m_verts[faceVerts[f][0]])[i];
570 ASSERTL0(
false,
"Could not find matching vertex for the face");
582 elementAaxis_length += pow(elementAaxis[i], 2);
583 elementBaxis_length += pow(elementBaxis[i], 2);
584 faceAaxis_length += pow(faceAaxis[i], 2);
585 faceBaxis_length += pow(faceBaxis[i], 2);
588 elementAaxis_length = sqrt(elementAaxis_length);
589 elementBaxis_length = sqrt(elementBaxis_length);
590 faceAaxis_length = sqrt(faceAaxis_length);
591 faceBaxis_length = sqrt(faceBaxis_length);
597 dotproduct1 += elementAaxis[i] * faceAaxis[i];
600 NekDouble norm = fabs(dotproduct1) / elementAaxis_length / faceAaxis_length;
610 if (dotproduct1 < 0.0)
618 dotproduct2 += elementBaxis[i] * faceBaxis[i];
621 norm = fabs(dotproduct2) / elementBaxis_length / faceBaxis_length;
625 "These vectors should be parallel");
629 if (dotproduct2 < 0.0)
644 dotproduct1 += elementAaxis[i] * faceBaxis[i];
647 norm = fabs(dotproduct1) / elementAaxis_length / faceBaxis_length;
651 "These vectors should be parallel");
655 if (dotproduct1 < 0.0)
664 dotproduct2 += elementBaxis[i] * faceAaxis[i];
667 norm = fabs(dotproduct2) / elementBaxis_length / faceAaxis_length;
671 "These vectors should be parallel");
673 if (dotproduct2 < 0.0)
679 orientation = orientation + 5;
690 for (
int i = 0; i < 4; ++i)
692 m_faces[i]->Reset(curvedEdges, curvedFaces);
700 for (
int i = 0; i < 4; ++i)
729 if (
m_xmap->GetBasisNumModes(0) != 2 ||
730 m_xmap->GetBasisNumModes(1) != 2 ||
731 m_xmap->GetBasisNumModes(2) != 2)
751 int order0 = *max_element(tmp.begin(), tmp.end());
754 tmp.push_back(order0);
757 int order1 = *max_element(tmp.begin(), tmp.end());
760 tmp.push_back(order0);
761 tmp.push_back(order1);
765 int order2 = *max_element(tmp.begin(), tmp.end());
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
#define ASSERTL0(condition, msg)
void SetUpFaceOrientation()
std::vector< StdRegions::Orientation > m_eorient
int GetEid(int i) const
Get the ID of edge i of this object.
Principle Modified Functions .
GeomFactorsSharedPtr m_geomFactors
Geometric factors.
StdRegions::StdExpansionSharedPtr GetXmap() const
Return the mapping object Geometry::m_xmap that represents the coordinate transformation from standar...
std::unordered_map< int, CurveSharedPtr > CurveMap
bool MinMaxCheck(const Array< OneD, const NekDouble > &gloCoord)
Check if given global coord is within twice the min/max distance of the element.
NekDouble dot(PointGeom &a)
retun the dot product between this and input a
int Imin(int n, const T *x, const int incx)
Return the index of the minimum element in x.
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 .
void Sub(PointGeom &a, PointGeom &b)
virtual void v_GenGeomFactors()
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
void ClampLocCoords(Array< OneD, NekDouble > &locCoord, NekDouble tol)
Clamp local coords to be within standard regions [-1, 1]^dim.
static const int kNtfaces
GeomState m_geomFactorsState
State of the geometric factors.
int GetGlobalID(void) const
Get the ID of this object.
std::vector< StdRegions::Orientation > m_forient
Gauss Radau pinned at x=-1, .
virtual int v_GetVertexFaceMap(const int i, const int j) const
Returns the standard element face IDs that are connected to a given vertex.
void SetUpEdgeOrientation()
static const unsigned int EdgeFaceConnectivity[6][2]
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)
std::shared_ptr< TriGeom > TriGeomSharedPtr
static const NekDouble kNekZeroTol
Principle Modified Functions .
static const unsigned int VertexFaceConnectivity[4][3]
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
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.
bool m_setupState
Wether or not the setup routines have been run.
int GetVid(int i) const
Get the ID of vertex i of this object.
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
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)
Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this ge...
void Mult(PointGeom &a, PointGeom &b)
_this = a x b
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
Geometry is straight-sided with constant geometric factors.
Gauss Radau pinned at x=-1, .
NekDouble GetLocCoords(const Array< OneD, const NekDouble > &coords, Array< OneD, NekDouble > &Lcoords)
Determine the local collapsed coordinates that correspond to a given Cartesian coordinate for this ge...
Geometry1DSharedPtr GetEdge(int i) const
Returns edge i of this object.
LibUtilities::ShapeType m_shapeType
Type of shape.
NekDouble dist(PointGeom &a)
return distance between this and input a
Geometric information has been generated.
void SetUpCoeffs(const int nCoeffs)
Initialise the Geometry::m_coeffs array.
static const unsigned int VertexEdgeConnectivity[4][3]
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
GeomType
Indicates the type of element geometry.
void SetUpLocalVertices()
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
std::shared_ptr< SegGeom > SegGeomSharedPtr
Geometry is curved or has non-constant factors.
PointGeomSharedPtr GetVertex(int i) const
Returns vertex i of this object.
int m_coordim
Coordinate dimension of this geometry object.
Describes the specification for a Basis.
virtual int v_GetVertexEdgeMap(const int i, const int j) const
Returns the standard element edge IDs that are connected to a given vertex.
virtual int v_GetEdgeFaceMap(const int i, const int j) const
Returns the standard element edge IDs that are connected to a given face.
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.
virtual int v_GetDir(const int faceidx, const int facedir) const
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &resid)
Determines if a point specified in global coordinates is located within this tetrahedral geometry and...