46 namespace SpatialDomains
49 const unsigned int HexGeom::VertexEdgeConnectivity[8][3] = {
50 {0, 3, 4}, {0, 1, 5}, {1, 2, 6}, {2, 3, 7},
51 {4, 8, 11}, {5, 8, 9}, {6, 9, 10}, {7, 10, 11}};
52 const unsigned int HexGeom::VertexFaceConnectivity[8][3] = {
53 {0, 1, 4}, {0, 1, 2}, {0, 2, 3}, {0, 3, 4},
54 {1, 4, 5}, {1, 2, 5}, {2, 3, 5}, {3, 4, 5}};
55 const unsigned int HexGeom::EdgeFaceConnectivity[12][2] = {
56 {0, 1}, {0, 2}, {0, 3}, {0, 4}, {1, 4}, {1, 2},
57 {2, 3}, {3, 4}, {1, 5}, {2, 5}, {3, 5}, {4, 5}};
65 :
Geometry3D(faces[0]->GetEdge(0)->GetVertex(0)->GetCoordim())
104 if (
m_xmap->GetBasisNumModes(0) != 2 ||
105 m_xmap->GetBasisNumModes(1) != 2 ||
106 m_xmap->GetBasisNumModes(2) != 2)
128 if (fabs((*
m_verts[faceVerts[f][0]])(i) -
129 (*
m_verts[faceVerts[f][1]])(i) +
130 (*
m_verts[faceVerts[f][2]])(i) -
131 (*
m_verts[faceVerts[f][3]])(i)) >
176 nq0 =
m_xmap->GetNumPoints(0);
177 nq1 =
m_xmap->GetNumPoints(1);
178 nq2 =
m_xmap->GetNumPoints(2);
184 len0 += (pts[nq0 - 1] - pts[0]) * (pts[nq0 - 1] - pts[0]);
185 xi0 += (coords[i] - pts[0]) * (pts[nq0 - 1] - pts[0]);
189 len1 += (pts[nq0 * (nq1 - 1)] - pts[0]) *
190 (pts[nq0 * (nq1 - 1)] - pts[0]);
191 xi1 += (coords[i] - pts[0]) * (pts[nq0 * (nq1 - 1)] - pts[0]);
195 len2 += (pts[nq0 * nq1 * (nq2 - 1)] - pts[0]) *
196 (pts[nq0 * nq1 * (nq2 - 1)] - pts[0]);
197 xi2 += (coords[i] - pts[0]) * (pts[nq0 * nq1 * (nq2 - 1)] - pts[0]);
200 Lcoords[0] = 2 * xi0 / len0 - 1.0;
201 Lcoords[1] = 2 * xi1 / len1 - 1.0;
202 Lcoords[2] = 2 * xi2 / len2 - 1.0;
206 PointGeom r(m_coordim, 0, coords[0], coords[1], coords[2]);
207 for (
int i = 0; i < 8; ++i)
215 int npts =
m_xmap->GetTotPoints();
238 ptdist = sqrt(tmp1[min_i]);
241 int qa = za.num_elements(), qb = zb.num_elements();
242 Lcoords[2] = zc[min_i / (qa * qb)];
243 min_i = min_i % (qa * qb);
244 Lcoords[1] = zb[min_i / qa];
245 Lcoords[0] = za[min_i % qa];
259 boost::ignore_unused(resid);
274 if (locCoord[0] >= -(1 + tol) && locCoord[0] <= 1 + tol &&
275 locCoord[1] >= -(1 + tol) && locCoord[1] <= 1 + tol &&
276 locCoord[2] >= -(1 + tol) && locCoord[2] <= 1 + tol)
304 if (faceidx == 0 || faceidx == 5)
308 else if (faceidx == 1 || faceidx == 3)
328 for (f = 1; f < 5; f++)
331 for (i = 0; i < 4; i++)
333 for (j = 0; j < 4; j++)
337 edge = std::dynamic_pointer_cast<
SegGeom>(
347 std::ostringstream errstrm;
348 errstrm <<
"Connected faces do not share an edge. Faces ";
355 std::ostringstream errstrm;
356 errstrm <<
"Connected faces share more than one edge. Faces ";
365 for (i = 0; i < 4; i++)
367 for (j = 0; j < 4; j++)
371 edge = std::dynamic_pointer_cast<
SegGeom>(
380 std::ostringstream errstrm;
381 errstrm <<
"Connected faces do not share an edge. Faces ";
388 std::ostringstream errstrm;
389 errstrm <<
"Connected faces share more than one edge. Faces ";
394 for (f = 1; f < 4; f++)
397 for (i = 0; i < 4; i++)
399 for (j = 0; j < 4; j++)
403 edge = std::dynamic_pointer_cast<
SegGeom>(
413 std::ostringstream errstrm;
414 errstrm <<
"Connected faces do not share an edge. Faces ";
421 std::ostringstream errstrm;
422 errstrm <<
"Connected faces share more than one edge. Faces ";
430 for (f = 1; f < 5; f++)
433 for (i = 0; i < 4; i++)
435 for (j = 0; j < 4; j++)
439 edge = std::dynamic_pointer_cast<
SegGeom>(
449 std::ostringstream errstrm;
450 errstrm <<
"Connected faces do not share an edge. Faces ";
457 std::ostringstream errstrm;
458 errstrm <<
"Connected faces share more than one edge. Faces ";
483 std::ostringstream errstrm;
484 errstrm <<
"Connected edges do not share a vertex. Edges ";
485 errstrm <<
m_edges[0]->GetGlobalID() <<
", " 492 for (i = 1; i < 3; i++)
504 std::ostringstream errstrm;
505 errstrm <<
"Connected edges do not share a vertex. Edges ";
506 errstrm <<
m_edges[i]->GetGlobalID() <<
", " 507 <<
m_edges[i - 1]->GetGlobalID();
528 std::ostringstream errstrm;
529 errstrm <<
"Connected edges do not share a vertex. Edges ";
530 errstrm <<
m_edges[8]->GetGlobalID() <<
", " 536 for (i = 9; i < 11; i++)
548 std::ostringstream errstrm;
549 errstrm <<
"Connected edges do not share a vertex. Edges ";
550 errstrm <<
m_edges[i]->GetGlobalID() <<
", " 551 <<
m_edges[i - 1]->GetGlobalID();
580 unsigned int baseVertex;
602 unsigned int orientation;
608 elementAaxis_length = 0.0;
609 elementBaxis_length = 0.0;
610 faceAaxis_length = 0.0;
611 faceBaxis_length = 0.0;
616 baseVertex =
m_faces[f]->GetVid(0);
632 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
633 (*
m_verts[faceVerts[f][0]])[i];
634 elementBaxis[i] = (*
m_verts[faceVerts[f][3]])[i] -
635 (*
m_verts[faceVerts[f][0]])[i];
642 elementAaxis[i] = (*
m_verts[faceVerts[f][1]])[i] -
643 (*
m_verts[faceVerts[f][0]])[i];
644 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
645 (*
m_verts[faceVerts[f][1]])[i];
652 elementAaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
653 (*
m_verts[faceVerts[f][3]])[i];
654 elementBaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
655 (*
m_verts[faceVerts[f][1]])[i];
662 elementAaxis[i] = (*
m_verts[faceVerts[f][2]])[i] -
663 (*
m_verts[faceVerts[f][3]])[i];
664 elementBaxis[i] = (*
m_verts[faceVerts[f][3]])[i] -
665 (*
m_verts[faceVerts[f][0]])[i];
670 ASSERTL0(
false,
"Could not find matching vertex for the face");
682 elementAaxis_length += pow(elementAaxis[i], 2);
683 elementBaxis_length += pow(elementBaxis[i], 2);
684 faceAaxis_length += pow(faceAaxis[i], 2);
685 faceBaxis_length += pow(faceBaxis[i], 2);
688 elementAaxis_length = sqrt(elementAaxis_length);
689 elementBaxis_length = sqrt(elementBaxis_length);
690 faceAaxis_length = sqrt(faceAaxis_length);
691 faceBaxis_length = sqrt(faceBaxis_length);
697 dotproduct1 += elementAaxis[i] * faceAaxis[i];
700 NekDouble norm = fabs(dotproduct1) / elementAaxis_length / faceAaxis_length;
710 if (dotproduct1 < 0.0)
718 dotproduct2 += elementBaxis[i] * faceBaxis[i];
721 norm = fabs(dotproduct2) / elementBaxis_length / faceBaxis_length;
725 "These vectors should be parallel");
729 if (dotproduct2 < 0.0)
744 dotproduct1 += elementAaxis[i] * faceBaxis[i];
747 norm = fabs(dotproduct1) / elementAaxis_length / faceBaxis_length;
751 "These vectors should be parallel");
755 if (dotproduct1 < 0.0)
764 dotproduct2 += elementBaxis[i] * faceAaxis[i];
767 norm = fabs(dotproduct2) / elementBaxis_length / faceAaxis_length;
771 "These vectors should be parallel");
773 if (dotproduct2 < 0.0)
779 orientation = orientation + 5;
791 const unsigned int edgeVerts[
kNedges][2] = {{0, 1},
817 ASSERTL0(
false,
"Could not find matching vertex for the edge");
826 for (
int i = 0; i < 6; ++i)
828 m_faces[i]->Reset(curvedEdges, curvedFaces);
839 for (
int i = 0; i < 6; ++i)
881 int order0 = *max_element(tmp1.begin(), tmp1.end());
907 int order1 = *max_element(tmp1.begin(), tmp1.end());
933 int order2 = *max_element(tmp1.begin(), tmp1.end());
static const unsigned int VertexEdgeConnectivity[8][3]
StdRegions::StdExpansionSharedPtr m_xmap
mapping containing isoparametric transformation.
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
int GetEid(int i) const
Get the ID of edge i of this object.
virtual void v_GenGeomFactors()
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
static const int kNtfaces
bool MinMaxCheck(const Array< OneD, const NekDouble > &gloCoord)
Check if given global coord is within twice the min/max distance of the element.
std::shared_ptr< QuadGeom > QuadGeomSharedPtr
virtual int v_GetEdgeFaceMap(const int i, const int j) const
Returns the standard element edge IDs that are connected to a given face.
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 .
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.
GeomState m_geomFactorsState
State of the geometric factors.
int GetGlobalID(void) const
Get the ID of this object.
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
void SetUpLocalVertices()
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...
virtual bool v_ContainsPoint(const Array< OneD, const NekDouble > &gloCoord, Array< OneD, NekDouble > &locCoord, NekDouble tol, NekDouble &resid)
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
virtual void v_Reset(CurveMap &curvedEdges, CurveMap &curvedFaces)
Reset this geometry object: unset the current state, zero Geometry::m_coeffs and remove allocated Geo...
virtual int v_GetVertexFaceMap(const int i, const int j) const
Returns the standard element face IDs that are connected to a given vertex.
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.
void SetUpEdgeOrientation()
virtual void v_FillGeom()
Put all quadrature information into face/edge structure and backward transform.
Array< OneD, Array< OneD, NekDouble > > m_coeffs
Array containing expansion coefficients of m_xmap.
static const unsigned int EdgeFaceConnectivity[12][2]
Geometry is straight-sided with constant geometric factors.
void SetUpFaceOrientation()
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.
static const unsigned int VertexFaceConnectivity[8][3]
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 int kNqfaces
GeomFactorsSharedPtr GetMetricInfo()
Get the geometric factors for this object.
virtual int v_GetDir(const int faceidx, const int facedir) const
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...
std::shared_ptr< SegGeom > SegGeomSharedPtr
Geometry is curved or has non-constant factors.
virtual int v_GetVertexEdgeMap(const int i, const int j) const
Returns the standard element edge IDs that are connected to a given vertex.
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.
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.